957a2454b45f76118c4d090f783001da084a240f
[alexxy/gromacs.git] / src / gromacs / gmxpreprocess / gen_ad.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
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.
10  *
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.
15  *
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.
20  *
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.
25  *
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.
33  *
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.
36  */
37 /* This file is completely threadsafe - keep it that way! */
38 #include "gmxpre.h"
39
40 #include "gen_ad.h"
41
42 #include <ctype.h>
43 #include <stdlib.h>
44 #include <string.h>
45
46 #include <cmath>
47
48 #include <algorithm>
49
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"
62
63 #define DIHEDRAL_WAS_SET_IN_RTP 0
64 static bool was_dihedral_set_in_rtp(const t_param *dih)
65 {
66     return dih->c[MAXFORCEPARAM-1] == DIHEDRAL_WAS_SET_IN_RTP;
67 }
68
69 typedef bool (*peq)(t_param *p1, t_param *p2);
70
71 static int acomp(const void *a1, const void *a2)
72 {
73     const t_param *p1, *p2;
74     int            ac;
75
76     p1 = static_cast<const t_param *>(a1);
77     p2 = static_cast<const t_param *>(a2);
78     if ((ac = (p1->aj()-p2->aj())) != 0)
79     {
80         return ac;
81     }
82     else if ((ac = (p1->ai()-p2->ai())) != 0)
83     {
84         return ac;
85     }
86     else
87     {
88         return (p1->ak()-p2->ak());
89     }
90 }
91
92 static int pcomp(const void *a1, const void *a2)
93 {
94     const t_param *p1, *p2;
95     int            pc;
96
97     p1 = static_cast<const t_param *>(a1);
98     p2 = static_cast<const t_param *>(a2);
99     if ((pc = (p1->ai()-p2->ai())) != 0)
100     {
101         return pc;
102     }
103     else
104     {
105         return (p1->aj()-p2->aj());
106     }
107 }
108
109 static int dcomp(const void *d1, const void *d2)
110 {
111     const t_param *p1, *p2;
112     int            dc;
113
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)
118     {
119         return dc;
120     }
121     else if ((dc = (p1->ak()-p2->ak())) != 0)
122     {
123         return dc;
124     }
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))
128     {
129         return -1;
130     }
131     else if (!was_dihedral_set_in_rtp(p1) &&
132              was_dihedral_set_in_rtp(p2))
133     {
134         return 1;
135     }
136     /* Then sort by I and J (two outer) atoms */
137     else if ((dc = (p1->ai()-p2->ai())) != 0)
138     {
139         return dc;
140     }
141     else if ((dc = (p1->al()-p2->al())) != 0)
142     {
143         return dc;
144     }
145     else
146     {
147         // AMBER force fields with type 9 dihedrals can reach here, where we sort on
148         // the contents of the string that names the macro for the parameters.
149         return strcmp(p1->s, p2->s);
150     }
151 }
152
153
154 static bool is_dihedral_on_same_bond(t_param *p1, t_param *p2)
155 {
156     return ((p1->aj() == p2->aj()) && (p1->ak() == p2->ak())) ||
157            ((p1->aj() == p2->ak()) && (p1->ak() == p2->aj()));
158 }
159
160
161 static bool preq(t_param *p1, t_param *p2)
162 {
163     return (p1->ai() == p2->ai()) && (p1->aj() == p2->aj());
164 }
165
166 static void rm2par(t_param p[], int *np, peq eq)
167 {
168     int *index, nind;
169     int  i, j;
170
171     if ((*np) == 0)
172     {
173         return;
174     }
175
176     snew(index, *np);
177     nind          = 0;
178     index[nind++] = 0;
179     for (i = 1; (i < (*np)); i++)
180     {
181         if (!eq(&p[i], &p[i-1]))
182         {
183             index[nind++] = i;
184         }
185     }
186     /* Index now holds pointers to all the non-equal params,
187      * this only works when p is sorted of course
188      */
189     for (i = 0; (i < nind); i++)
190     {
191         for (j = 0; (j < MAXATOMLIST); j++)
192         {
193             p[i].a[j] = p[index[i]].a[j];
194         }
195         for (j = 0; (j < MAXFORCEPARAM); j++)
196         {
197             p[i].c[j] = p[index[i]].c[j];
198         }
199         if (p[index[i]].a[0] == p[index[i]].a[1])
200         {
201             strcpy(p[i].s, "");
202         }
203         else if (index[i] > i)
204         {
205             /* Copy the string only if it comes from somewhere else
206              * otherwise we will end up copying a random (newly freed) pointer.
207              * Since the index is sorted we only have to test for index[i] > i.
208              */
209             strcpy(p[i].s, p[index[i]].s);
210         }
211     }
212     (*np) = nind;
213
214     sfree(index);
215 }
216
217 static void cppar(t_param p[], int np, t_params plist[], int ftype)
218 {
219     int       i, j, nral, nrfp;
220     t_params *ps;
221
222     ps   = &plist[ftype];
223     nral = NRAL(ftype);
224     nrfp = NRFP(ftype);
225
226     /* Keep old stuff */
227     pr_alloc(np, ps);
228     for (i = 0; (i < np); i++)
229     {
230         for (j = 0; (j < nral); j++)
231         {
232             ps->param[ps->nr].a[j] = p[i].a[j];
233         }
234         for (j = 0; (j < nrfp); j++)
235         {
236             ps->param[ps->nr].c[j] = p[i].c[j];
237         }
238         for (j = 0; (j < MAXSLEN); j++)
239         {
240             ps->param[ps->nr].s[j] = p[i].s[j];
241         }
242         ps->nr++;
243     }
244 }
245
246 static void cpparam(t_param *dest, t_param *src)
247 {
248     int j;
249
250     for (j = 0; (j < MAXATOMLIST); j++)
251     {
252         dest->a[j] = src->a[j];
253     }
254     for (j = 0; (j < MAXFORCEPARAM); j++)
255     {
256         dest->c[j] = src->c[j];
257     }
258     for (j = 0; (j < MAXSLEN); j++)
259     {
260         dest->s[j] = src->s[j];
261     }
262 }
263
264 static void set_p(t_param *p, const int ai[4], const real *c, char *s)
265 {
266     int j;
267
268     for (j = 0; (j < 4); j++)
269     {
270         p->a[j] = ai[j];
271     }
272     for (j = 0; (j < MAXFORCEPARAM); j++)
273     {
274         if (c)
275         {
276             p->c[j] = c[j];
277         }
278         else
279         {
280             p->c[j] = NOTSET;
281         }
282     }
283
284     set_p_string(p, s);
285 }
286
287 static int idcomp(const void *a, const void *b)
288 {
289     const t_param *pa, *pb;
290     int            d;
291
292     pa = static_cast<const t_param *>(a);
293     pb = static_cast<const t_param *>(b);
294     if ((d = (pa->a[0]-pb->a[0])) != 0)
295     {
296         return d;
297     }
298     else if ((d = (pa->a[3]-pb->a[3])) != 0)
299     {
300         return d;
301     }
302     else if ((d = (pa->a[1]-pb->a[1])) != 0)
303     {
304         return d;
305     }
306     else
307     {
308         return (pa->a[2]-pb->a[2]);
309     }
310 }
311
312 static void sort_id(int nr, t_param ps[])
313 {
314     int i, tmp;
315
316     /* First swap order of atoms around if necessary */
317     for (i = 0; (i < nr); i++)
318     {
319         if (ps[i].a[3] < ps[i].a[0])
320         {
321             tmp = ps[i].a[3]; ps[i].a[3] = ps[i].a[0]; ps[i].a[0] = tmp;
322             tmp = ps[i].a[2]; ps[i].a[2] = ps[i].a[1]; ps[i].a[1] = tmp;
323         }
324     }
325     /* Now sort it */
326     if (nr > 1)
327     {
328         qsort(ps, nr, static_cast<size_t>(sizeof(ps[0])), idcomp);
329     }
330 }
331
332 static int n_hydro(const int a[], char ***atomname)
333 {
334     int  i, nh = 0;
335     char c0, c1, *aname;
336
337     for (i = 0; (i < 4); i += 3)
338     {
339         aname = *atomname[a[i]];
340         c0    = toupper(aname[0]);
341         if (c0 == 'H')
342         {
343             nh++;
344         }
345         else if ((static_cast<int>(strlen(aname)) > 1) && (c0 >= '0') && (c0 <= '9'))
346         {
347             c1 = toupper(aname[1]);
348             if (c1 == 'H')
349             {
350                 nh++;
351             }
352         }
353     }
354     return nh;
355 }
356
357 /* Clean up the dihedrals (both generated and read from the .rtp
358  * file). */
359 static void clean_dih(t_param *dih, int *ndih, t_param improper[], int nimproper,
360                       t_atoms *atoms, bool bKeepAllGeneratedDihedrals,
361                       bool bRemoveDihedralIfWithImproper)
362 {
363     int   i, j, k, l;
364     int  *index, nind;
365
366     /* Construct the list of the indices of the dihedrals
367      * (i.e. generated or read) that might be kept. */
368     snew(index, *ndih+1);
369     if (bKeepAllGeneratedDihedrals)
370     {
371         fprintf(stderr, "Keeping all generated dihedrals\n");
372         nind = *ndih;
373         for (i = 0; i < nind; i++)
374         {
375             index[i] = i;
376         }
377         index[nind] = *ndih;
378     }
379     else
380     {
381         nind = 0;
382         /* Check if generated dihedral i should be removed. The
383          * dihedrals have been sorted by dcomp() above, so all those
384          * on the same two central atoms are together, with those from
385          * the .rtp file preceding those that were automatically
386          * generated. We remove the latter if the former exist. */
387         for (i = 0; i < *ndih; i++)
388         {
389             /* Keep the dihedrals that were defined in the .rtp file,
390              * and the dihedrals that were generated and different
391              * from the last one (whether it was generated or not). */
392             if (was_dihedral_set_in_rtp(&dih[i]) ||
393                 0 == i ||
394                 !is_dihedral_on_same_bond(&dih[i], &dih[i-1]))
395             {
396                 index[nind++] = i;
397             }
398         }
399         index[nind] = *ndih;
400     }
401
402     k = 0;
403     for (i = 0; i < nind; i++)
404     {
405         bool bWasSetInRTP = was_dihedral_set_in_rtp(&dih[index[i]]);
406         bool bKeep        = TRUE;
407         if (!bWasSetInRTP && bRemoveDihedralIfWithImproper)
408         {
409             /* Remove the dihedral if there is an improper on the same
410              * bond. */
411             for (j = 0; j < nimproper && bKeep; j++)
412             {
413                 bKeep = !is_dihedral_on_same_bond(&dih[index[i]], &improper[j]);
414             }
415         }
416
417         if (bKeep)
418         {
419             /* If we don't want all dihedrals, we want to select the
420              * ones with the fewest hydrogens. Note that any generated
421              * dihedrals on the same bond as an .rtp dihedral may have
422              * been already pruned above in the construction of
423              * index[]. However, their parameters are still present,
424              * and l is looping over this dihedral and all of its
425              * pruned siblings. */
426             int bestl = index[i];
427             if (!bKeepAllGeneratedDihedrals && !bWasSetInRTP)
428             {
429                 /* Minimum number of hydrogens for i and l atoms */
430                 int minh = 2;
431                 for (l = index[i];
432                      (l < index[i+1] &&
433                       is_dihedral_on_same_bond(&dih[index[i]], &dih[l]));
434                      l++)
435                 {
436                     int nh = n_hydro(dih[l].a, atoms->atomname);
437                     if (nh < minh)
438                     {
439                         minh  = nh;
440                         bestl = l;
441                     }
442                     if (0 == minh)
443                     {
444                         break;
445                     }
446                 }
447             }
448             if (k != bestl)
449             {
450                 cpparam(&dih[k], &dih[bestl]);
451             }
452             k++;
453         }
454     }
455
456     for (i = k; i < *ndih; i++)
457     {
458         strcpy(dih[i].s, "");
459     }
460     *ndih = k;
461
462     sfree(index);
463 }
464
465 static int get_impropers(t_atoms *atoms, t_hackblock hb[], t_param **improper,
466                          bool bAllowMissing)
467 {
468     t_rbondeds   *impropers;
469     int           nimproper, i, j, k, start, ninc, nalloc;
470     int           ai[MAXATOMLIST];
471     bool          bStop;
472
473     ninc   = 500;
474     nalloc = ninc;
475     snew(*improper, nalloc);
476
477     /* Add all the impropers from the residue database to the list. */
478     nimproper = 0;
479     start     = 0;
480     if (hb != nullptr)
481     {
482         for (i = 0; (i < atoms->nres); i++)
483         {
484             impropers = &hb[i].rb[ebtsIDIHS];
485             for (j = 0; (j < impropers->nb); j++)
486             {
487                 bStop = FALSE;
488                 for (k = 0; (k < 4) && !bStop; k++)
489                 {
490                     ai[k] = search_atom(impropers->b[j].a[k], start,
491                                         atoms,
492                                         "improper", bAllowMissing);
493                     if (ai[k] == -1)
494                     {
495                         bStop = TRUE;
496                     }
497                 }
498                 if (!bStop)
499                 {
500                     if (nimproper == nalloc)
501                     {
502                         nalloc += ninc;
503                         srenew(*improper, nalloc);
504                     }
505                     /* Not broken out */
506                     set_p(&((*improper)[nimproper]), ai, nullptr, impropers->b[j].s);
507                     nimproper++;
508                 }
509             }
510             while ((start < atoms->nr) && (atoms->atom[start].resind == i))
511             {
512                 start++;
513             }
514         }
515     }
516
517     return nimproper;
518 }
519
520 static int nb_dist(t_nextnb *nnb, int ai, int aj)
521 {
522     int  nre, nrx, NRE;
523     int *nrexcl;
524     int *a;
525
526     if (ai == aj)
527     {
528         return 0;
529     }
530
531     NRE    = -1;
532     nrexcl = nnb->nrexcl[ai];
533     for (nre = 1; (nre < nnb->nrex); nre++)
534     {
535         a = nnb->a[ai][nre];
536         for (nrx = 0; (nrx < nrexcl[nre]); nrx++)
537         {
538             if ((aj == a[nrx]) && (NRE == -1))
539             {
540                 NRE = nre;
541             }
542         }
543     }
544     return NRE;
545 }
546
547 static bool is_hydro(t_atoms *atoms, int ai)
548 {
549     return ((*(atoms->atomname[ai]))[0] == 'H');
550 }
551
552 static void get_atomnames_min(int n, char **anm,
553                               int resind, t_atoms *atoms, const int *a)
554 {
555     int m;
556
557     /* Assume ascending residue numbering */
558     for (m = 0; m < n; m++)
559     {
560         if (atoms->atom[a[m]].resind < resind)
561         {
562             strcpy(anm[m], "-");
563         }
564         else if (atoms->atom[a[m]].resind > resind)
565         {
566             strcpy(anm[m], "+");
567         }
568         else
569         {
570             strcpy(anm[m], "");
571         }
572         strcat(anm[m], *(atoms->atomname[a[m]]));
573     }
574 }
575
576 static void gen_excls(t_atoms *atoms, t_excls *excls, t_hackblock hb[],
577                       bool bAllowMissing)
578 {
579     int         r;
580     int         a, astart, i1, i2, itmp;
581     t_rbondeds *hbexcl;
582     int         e;
583     char       *anm;
584
585     astart = 0;
586     for (a = 0; a < atoms->nr; a++)
587     {
588         r = atoms->atom[a].resind;
589         if (a == atoms->nr-1 || atoms->atom[a+1].resind != r)
590         {
591             hbexcl = &hb[r].rb[ebtsEXCLS];
592
593             for (e = 0; e < hbexcl->nb; e++)
594             {
595                 anm = hbexcl->b[e].a[0];
596                 i1  = search_atom(anm, astart, atoms,
597                                   "exclusion", bAllowMissing);
598                 anm = hbexcl->b[e].a[1];
599                 i2  = search_atom(anm, astart, atoms,
600                                   "exclusion", bAllowMissing);
601                 if (i1 != -1 && i2 != -1)
602                 {
603                     if (i1 > i2)
604                     {
605                         itmp = i1;
606                         i1   = i2;
607                         i2   = itmp;
608                     }
609                     srenew(excls[i1].e, excls[i1].nr+1);
610                     excls[i1].e[excls[i1].nr] = i2;
611                     excls[i1].nr++;
612                 }
613             }
614
615             astart = a+1;
616         }
617     }
618
619     for (a = 0; a < atoms->nr; a++)
620     {
621         if (excls[a].nr > 1)
622         {
623             std::sort(excls[a].e, excls[a].e+excls[a].nr);
624         }
625     }
626 }
627
628 static void remove_excl(t_excls *excls, int remove)
629 {
630     int i;
631
632     for (i = remove+1; i < excls->nr; i++)
633     {
634         excls->e[i-1] = excls->e[i];
635     }
636
637     excls->nr--;
638 }
639
640 void clean_excls(t_nextnb *nnb, int nrexcl, t_excls excls[])
641 {
642     int      i, j, j1, k, k1, l, l1, e;
643     t_excls *excl;
644
645     if (nrexcl >= 1)
646     {
647         /* extract all i-j-k-l neighbours from nnb struct */
648         for (i = 0; (i < nnb->nr); i++)
649         {
650             /* For all particles */
651             excl = &excls[i];
652
653             for (j = 0; (j < nnb->nrexcl[i][1]); j++)
654             {
655                 /* For all first neighbours */
656                 j1 = nnb->a[i][1][j];
657
658                 for (e = 0; e < excl->nr; e++)
659                 {
660                     if (excl->e[e] == j1)
661                     {
662                         remove_excl(excl, e);
663                     }
664                 }
665
666                 if (nrexcl >= 2)
667                 {
668                     for (k = 0; (k < nnb->nrexcl[j1][1]); k++)
669                     {
670                         /* For all first neighbours of j1 */
671                         k1 = nnb->a[j1][1][k];
672
673                         for (e = 0; e < excl->nr; e++)
674                         {
675                             if (excl->e[e] == k1)
676                             {
677                                 remove_excl(excl, e);
678                             }
679                         }
680
681                         if (nrexcl >= 3)
682                         {
683                             for (l = 0; (l < nnb->nrexcl[k1][1]); l++)
684                             {
685                                 /* For all first neighbours of k1 */
686                                 l1 = nnb->a[k1][1][l];
687
688                                 for (e = 0; e < excl->nr; e++)
689                                 {
690                                     if (excl->e[e] == l1)
691                                     {
692                                         remove_excl(excl, e);
693                                     }
694                                 }
695                             }
696                         }
697                     }
698                 }
699             }
700         }
701     }
702 }
703
704 void generate_excls(t_nextnb *nnb, int nrexcl, t_excls excls[])
705 {
706     int      i, j, n, N;
707     t_excls *excl;
708
709     for (N = 1; (N < std::min(nrexcl, nnb->nrex)); N++)
710     {
711         /* extract all i-j-k-l neighbours from nnb struct */
712         for (i = 0; (i < nnb->nr); i++)
713         {
714             /* For all particles */
715             excl      = &excls[i];
716             n         = excl->nr;
717             excl->nr += nnb->nrexcl[i][N];
718             srenew(excl->e, excl->nr);
719             for (j = 0; (j < nnb->nrexcl[i][N]); j++)
720             {
721                 /* For all first neighbours */
722                 if (nnb->a[i][N][j] != i)
723                 {
724                     excl->e[n++] = nnb->a[i][N][j];
725                 }
726             }
727         }
728     }
729 }
730
731 /* Generate pairs, angles and dihedrals from .rtp settings */
732 void gen_pad(t_nextnb *nnb, t_atoms *atoms, t_restp rtp[],
733              t_params plist[], t_excls excls[], t_hackblock hb[],
734              bool bAllowMissing)
735 {
736     t_param    *ang, *dih, *pai, *improper;
737     t_rbondeds *hbang, *hbdih;
738     char      **anm;
739     const char *p;
740     int         res, minres, maxres;
741     int         i, j, j1, k, k1, l, l1, m, n, i1, i2;
742     int         ninc, maxang, maxdih, maxpai;
743     int         nang, ndih, npai, nimproper, nbd;
744     int         nFound;
745     bool        bFound, bExcl;
746
747     /* These are the angles, dihedrals and pairs that we generate
748      * from the bonds. The ones that are already there from the rtp file
749      * will be retained.
750      */
751     nang   = 0;
752     npai   = 0;
753     ndih   = 0;
754     ninc   = 500;
755     maxang = maxdih = maxpai = ninc;
756     snew(ang, maxang);
757     snew(dih, maxdih);
758     snew(pai, maxpai);
759
760     snew(anm, 4);
761     for (i = 0; i < 4; i++)
762     {
763         snew(anm[i], 12);
764     }
765
766     if (hb)
767     {
768         gen_excls(atoms, excls, hb, bAllowMissing);
769         /* mark all entries as not matched yet */
770         for (i = 0; i < atoms->nres; i++)
771         {
772             for (j = 0; j < ebtsNR; j++)
773             {
774                 for (k = 0; k < hb[i].rb[j].nb; k++)
775                 {
776                     hb[i].rb[j].b[k].match = FALSE;
777                 }
778             }
779         }
780     }
781
782     /* Extract all i-j-k-l neighbours from nnb struct to generate all
783      * angles and dihedrals. */
784     for (i = 0; (i < nnb->nr); i++)
785     {
786         /* For all particles */
787         for (j = 0; (j < nnb->nrexcl[i][1]); j++)
788         {
789             /* For all first neighbours */
790             j1 = nnb->a[i][1][j];
791             for (k = 0; (k < nnb->nrexcl[j1][1]); k++)
792             {
793                 /* For all first neighbours of j1 */
794                 k1 = nnb->a[j1][1][k];
795                 if (k1 != i)
796                 {
797                     /* Generate every angle only once */
798                     if (i < k1)
799                     {
800                         if (nang == maxang)
801                         {
802                             maxang += ninc;
803                             srenew(ang, maxang);
804                         }
805                         ang[nang].ai() = i;
806                         ang[nang].aj() = j1;
807                         ang[nang].ak() = k1;
808                         ang[nang].c0() = NOTSET;
809                         ang[nang].c1() = NOTSET;
810                         set_p_string(&(ang[nang]), "");
811                         if (hb)
812                         {
813                             minres = atoms->atom[ang[nang].a[0]].resind;
814                             maxres = minres;
815                             for (m = 1; m < 3; m++)
816                             {
817                                 minres = std::min(minres, atoms->atom[ang[nang].a[m]].resind);
818                                 maxres = std::max(maxres, atoms->atom[ang[nang].a[m]].resind);
819                             }
820                             res = 2*minres-maxres;
821                             do
822                             {
823                                 res += maxres-minres;
824                                 get_atomnames_min(3, anm, res, atoms, ang[nang].a);
825                                 hbang = &hb[res].rb[ebtsANGLES];
826                                 for (l = 0; (l < hbang->nb); l++)
827                                 {
828                                     if (strcmp(anm[1], hbang->b[l].aj()) == 0)
829                                     {
830                                         bFound = FALSE;
831                                         for (m = 0; m < 3; m += 2)
832                                         {
833                                             bFound = (bFound ||
834                                                       ((strcmp(anm[m], hbang->b[l].ai()) == 0) &&
835                                                        (strcmp(anm[2-m], hbang->b[l].ak()) == 0)));
836                                         }
837                                         if (bFound)
838                                         {
839                                             set_p_string(&(ang[nang]), hbang->b[l].s);
840                                             /* Mark that we found a match for this entry */
841                                             hbang->b[l].match = TRUE;
842                                         }
843                                     }
844                                 }
845                             }
846                             while (res < maxres);
847                         }
848                         nang++;
849                     }
850                     /* Generate every dihedral, 1-4 exclusion and 1-4 interaction
851                        only once */
852                     if (j1 < k1)
853                     {
854                         for (l = 0; (l < nnb->nrexcl[k1][1]); l++)
855                         {
856                             /* For all first neighbours of k1 */
857                             l1 = nnb->a[k1][1][l];
858                             if ((l1 != i) && (l1 != j1))
859                             {
860                                 if (ndih == maxdih)
861                                 {
862                                     maxdih += ninc;
863                                     srenew(dih, maxdih);
864                                 }
865                                 dih[ndih].ai() = i;
866                                 dih[ndih].aj() = j1;
867                                 dih[ndih].ak() = k1;
868                                 dih[ndih].al() = l1;
869                                 for (m = 0; m < MAXFORCEPARAM; m++)
870                                 {
871                                     dih[ndih].c[m] = NOTSET;
872                                 }
873                                 set_p_string(&(dih[ndih]), "");
874                                 nFound = 0;
875                                 if (hb)
876                                 {
877                                     minres = atoms->atom[dih[ndih].a[0]].resind;
878                                     maxres = minres;
879                                     for (m = 1; m < 4; m++)
880                                     {
881                                         minres = std::min(minres, atoms->atom[dih[ndih].a[m]].resind);
882                                         maxres = std::max(maxres, atoms->atom[dih[ndih].a[m]].resind);
883                                     }
884                                     res = 2*minres-maxres;
885                                     do
886                                     {
887                                         res += maxres-minres;
888                                         get_atomnames_min(4, anm, res, atoms, dih[ndih].a);
889                                         hbdih = &hb[res].rb[ebtsPDIHS];
890                                         for (n = 0; (n < hbdih->nb); n++)
891                                         {
892                                             bFound = FALSE;
893                                             for (m = 0; m < 2; m++)
894                                             {
895                                                 bFound = (bFound ||
896                                                           ((strcmp(anm[3*m],  hbdih->b[n].ai()) == 0) &&
897                                                            (strcmp(anm[1+m],  hbdih->b[n].aj()) == 0) &&
898                                                            (strcmp(anm[2-m],  hbdih->b[n].ak()) == 0) &&
899                                                            (strcmp(anm[3-3*m], hbdih->b[n].al()) == 0)));
900                                             }
901                                             if (bFound)
902                                             {
903                                                 set_p_string(&dih[ndih], hbdih->b[n].s);
904                                                 /* Mark that we found a match for this entry */
905                                                 hbdih->b[n].match = TRUE;
906
907                                                 /* Set the last parameter to be able to see
908                                                    if the dihedral was in the rtp list.
909                                                  */
910                                                 dih[ndih].c[MAXFORCEPARAM-1] = DIHEDRAL_WAS_SET_IN_RTP;
911                                                 nFound++;
912                                                 ndih++;
913                                                 /* Set the next direct in case the rtp contains
914                                                    multiple entries for this dihedral.
915                                                  */
916                                                 if (ndih == maxdih)
917                                                 {
918                                                     maxdih += ninc;
919                                                     srenew(dih, maxdih);
920                                                 }
921                                                 dih[ndih].ai() = i;
922                                                 dih[ndih].aj() = j1;
923                                                 dih[ndih].ak() = k1;
924                                                 dih[ndih].al() = l1;
925                                                 for (m = 0; m < MAXFORCEPARAM; m++)
926                                                 {
927                                                     dih[ndih].c[m] = NOTSET;
928                                                 }
929                                             }
930                                         }
931                                     }
932                                     while (res < maxres);
933                                 }
934                                 if (nFound == 0)
935                                 {
936                                     if (ndih == maxdih)
937                                     {
938                                         maxdih += ninc;
939                                         srenew(dih, maxdih);
940                                     }
941                                     dih[ndih].ai() = i;
942                                     dih[ndih].aj() = j1;
943                                     dih[ndih].ak() = k1;
944                                     dih[ndih].al() = l1;
945                                     for (m = 0; m < MAXFORCEPARAM; m++)
946                                     {
947                                         dih[ndih].c[m] = NOTSET;
948                                     }
949                                     set_p_string(&(dih[ndih]), "");
950                                     ndih++;
951                                 }
952
953                                 nbd = nb_dist(nnb, i, l1);
954                                 if (nbd == 3)
955                                 {
956                                     i1    = std::min(i, l1);
957                                     i2    = std::max(i, l1);
958                                     bExcl = FALSE;
959                                     for (m = 0; m < excls[i1].nr; m++)
960                                     {
961                                         bExcl = bExcl || excls[i1].e[m] == i2;
962                                     }
963                                     if (!bExcl)
964                                     {
965                                         if (rtp[0].bGenerateHH14Interactions ||
966                                             !(is_hydro(atoms, i1) && is_hydro(atoms, i2)))
967                                         {
968                                             if (npai == maxpai)
969                                             {
970                                                 maxpai += ninc;
971                                                 srenew(pai, maxpai);
972                                             }
973                                             pai[npai].ai() = i1;
974                                             pai[npai].aj() = i2;
975                                             pai[npai].c0() = NOTSET;
976                                             pai[npai].c1() = NOTSET;
977                                             set_p_string(&(pai[npai]), "");
978                                             npai++;
979                                         }
980                                     }
981                                 }
982                             }
983                         }
984                     }
985                 }
986             }
987         }
988     }
989
990     if (hb)
991     {
992         /* The above approach is great in that we double-check that e.g. an angle
993          * really corresponds to three atoms connected by bonds, but this is not
994          * generally true. Go through the angle and dihedral hackblocks to add
995          * entries that we have not yet marked as matched when going through bonds.
996          */
997         for (i = 0; i < atoms->nres; i++)
998         {
999             /* Add remaining angles from hackblock */
1000             hbang = &hb[i].rb[ebtsANGLES];
1001             for (j = 0; j < hbang->nb; j++)
1002             {
1003                 if (hbang->b[j].match == TRUE)
1004                 {
1005                     /* We already used this entry, continue to the next */
1006                     continue;
1007                 }
1008                 /* Hm - entry not used, let's see if we can find all atoms */
1009                 if (nang == maxang)
1010                 {
1011                     maxang += ninc;
1012                     srenew(ang, maxang);
1013                 }
1014                 bFound = TRUE;
1015                 for (k = 0; k < 3 && bFound; k++)
1016                 {
1017                     p   = hbang->b[j].a[k];
1018                     res = i;
1019                     if (p[0] == '-')
1020                     {
1021                         p++;
1022                         res--;
1023                     }
1024                     else if (p[0] == '+')
1025                     {
1026                         p++;
1027                         res++;
1028                     }
1029                     ang[nang].a[k] = search_res_atom(p, res, atoms, "angle", TRUE);
1030                     bFound         = (ang[nang].a[k] != -1);
1031                 }
1032                 ang[nang].c0() = NOTSET;
1033                 ang[nang].c1() = NOTSET;
1034
1035                 if (bFound)
1036                 {
1037                     set_p_string(&(ang[nang]), hbang->b[j].s);
1038                     hbang->b[j].match = TRUE;
1039                     /* Incrementing nang means we save this angle */
1040                     nang++;
1041                 }
1042             }
1043
1044             /* Add remaining dihedrals from hackblock */
1045             hbdih = &hb[i].rb[ebtsPDIHS];
1046             for (j = 0; j < hbdih->nb; j++)
1047             {
1048                 if (hbdih->b[j].match == TRUE)
1049                 {
1050                     /* We already used this entry, continue to the next */
1051                     continue;
1052                 }
1053                 /* Hm - entry not used, let's see if we can find all atoms */
1054                 if (ndih == maxdih)
1055                 {
1056                     maxdih += ninc;
1057                     srenew(dih, maxdih);
1058                 }
1059                 bFound = TRUE;
1060                 for (k = 0; k < 4 && bFound; k++)
1061                 {
1062                     p   = hbdih->b[j].a[k];
1063                     res = i;
1064                     if (p[0] == '-')
1065                     {
1066                         p++;
1067                         res--;
1068                     }
1069                     else if (p[0] == '+')
1070                     {
1071                         p++;
1072                         res++;
1073                     }
1074                     dih[ndih].a[k] = search_res_atom(p, res, atoms, "dihedral", TRUE);
1075                     bFound         = (dih[ndih].a[k] != -1);
1076                 }
1077                 for (m = 0; m < MAXFORCEPARAM; m++)
1078                 {
1079                     dih[ndih].c[m] = NOTSET;
1080                 }
1081
1082                 if (bFound)
1083                 {
1084                     set_p_string(&(dih[ndih]), hbdih->b[j].s);
1085                     hbdih->b[j].match = TRUE;
1086                     /* Incrementing ndih means we save this dihedral */
1087                     ndih++;
1088                 }
1089             }
1090         }
1091     }
1092
1093     /* Sort angles with respect to j-i-k (middle atom first) */
1094     if (nang > 1)
1095     {
1096         qsort(ang, nang, static_cast<size_t>(sizeof(ang[0])), acomp);
1097     }
1098
1099     /* Sort dihedrals with respect to j-k-i-l (middle atoms first) */
1100     if (ndih > 1)
1101     {
1102         qsort(dih, ndih, static_cast<size_t>(sizeof(dih[0])), dcomp);
1103     }
1104
1105     /* Sort the pairs */
1106     if (npai > 1)
1107     {
1108         qsort(pai, npai, static_cast<size_t>(sizeof(pai[0])), pcomp);
1109     }
1110     if (npai > 0)
1111     {
1112         /* Remove doubles, could occur in 6-rings, such as phenyls,
1113            maybe one does not want this when fudgeQQ < 1.
1114          */
1115         fprintf(stderr, "Before cleaning: %d pairs\n", npai);
1116         rm2par(pai, &npai, preq);
1117     }
1118
1119     /* Get the impropers from the database */
1120     nimproper = get_impropers(atoms, hb, &improper, bAllowMissing);
1121
1122     /* Sort the impropers */
1123     sort_id(nimproper, improper);
1124
1125     if (ndih > 0)
1126     {
1127         fprintf(stderr, "Before cleaning: %d dihedrals\n", ndih);
1128         clean_dih(dih, &ndih, improper, nimproper, atoms,
1129                   rtp[0].bKeepAllGeneratedDihedrals,
1130                   rtp[0].bRemoveDihedralIfWithImproper);
1131     }
1132
1133     /* Now we have unique lists of angles and dihedrals
1134      * Copy them into the destination struct
1135      */
1136     cppar(ang, nang, plist, F_ANGLES);
1137     cppar(dih, ndih, plist, F_PDIHS);
1138     cppar(improper, nimproper, plist, F_IDIHS);
1139     cppar(pai, npai, plist, F_LJ14);
1140
1141     /* Remove all exclusions which are within nrexcl */
1142     clean_excls(nnb, rtp[0].nrexcl, excls);
1143
1144     sfree(ang);
1145     sfree(dih);
1146     sfree(improper);
1147     sfree(pai);
1148 }