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