Merge release-5-0 into master
[alexxy/gromacs.git] / src / gromacs / topology / topsort.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-2008, 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 #include "gmxpre.h"
38
39 #include "gromacs/topology/topsort.h"
40
41 #include "config.h"
42
43 #include <stdio.h>
44
45 #include "gromacs/legacyheaders/types/ifunc.h"
46
47 #include "gromacs/topology/topology.h"
48 #include "gromacs/utility/fatalerror.h"
49 #include "gromacs/utility/smalloc.h"
50
51 static gmx_bool ip_pert(int ftype, const t_iparams *ip)
52 {
53     gmx_bool bPert;
54     int      i;
55
56     if (NRFPB(ftype) == 0)
57     {
58         return FALSE;
59     }
60
61     switch (ftype)
62     {
63         case F_BONDS:
64         case F_G96BONDS:
65         case F_HARMONIC:
66         case F_ANGLES:
67         case F_G96ANGLES:
68         case F_IDIHS:
69             bPert = (ip->harmonic.rA  != ip->harmonic.rB ||
70                      ip->harmonic.krA != ip->harmonic.krB);
71             break;
72         case F_MORSE:
73             bPert = (ip->morse.b0A  != ip->morse.b0B ||
74                      ip->morse.cbA  != ip->morse.cbB ||
75                      ip->morse.betaA  != ip->morse.betaB);
76             break;
77         case F_RESTRBONDS:
78             bPert = (ip->restraint.lowA  != ip->restraint.lowB ||
79                      ip->restraint.up1A  != ip->restraint.up1B ||
80                      ip->restraint.up2A  != ip->restraint.up2B ||
81                      ip->restraint.kA    != ip->restraint.kB);
82             break;
83         case F_UREY_BRADLEY:
84             bPert = (ip->u_b.thetaA  != ip->u_b.thetaB  ||
85                      ip->u_b.kthetaA != ip->u_b.kthetaB ||
86                      ip->u_b.r13A    != ip->u_b.r13B    ||
87                      ip->u_b.kUBA    != ip->u_b.kUBB);
88             break;
89         case F_PDIHS:
90         case F_PIDIHS:
91         case F_ANGRES:
92         case F_ANGRESZ:
93             bPert = (ip->pdihs.phiA != ip->pdihs.phiB ||
94                      ip->pdihs.cpA  != ip->pdihs.cpB);
95             break;
96         case F_RBDIHS:
97             bPert = FALSE;
98             for (i = 0; i < NR_RBDIHS; i++)
99             {
100                 if (ip->rbdihs.rbcA[i] != ip->rbdihs.rbcB[i])
101                 {
102                     bPert = TRUE;
103                 }
104             }
105             break;
106         case F_TABBONDS:
107         case F_TABBONDSNC:
108         case F_TABANGLES:
109         case F_TABDIHS:
110             bPert = (ip->tab.kA != ip->tab.kB);
111             break;
112         case F_POSRES:
113             bPert = FALSE;
114             for (i = 0; i < DIM; i++)
115             {
116                 if (ip->posres.pos0A[i] != ip->posres.pos0B[i] ||
117                     ip->posres.fcA[i]   != ip->posres.fcB[i])
118                 {
119                     bPert = TRUE;
120                 }
121             }
122             break;
123         case F_DIHRES:
124             bPert = ((ip->dihres.phiA != ip->dihres.phiB) ||
125                      (ip->dihres.dphiA != ip->dihres.dphiB) ||
126                      (ip->dihres.kfacA != ip->dihres.kfacB));
127             break;
128         case F_LJ14:
129             bPert = (ip->lj14.c6A  != ip->lj14.c6B ||
130                      ip->lj14.c12A != ip->lj14.c12B);
131             break;
132         case F_CMAP:
133             bPert = FALSE;
134             break;
135         case F_RESTRANGLES:
136         case F_RESTRDIHS:
137         case F_CBTDIHS:
138             bPert = FALSE;
139             gmx_fatal(FARGS, "Function type %s does not support currentely free energy calculations",
140                       interaction_function[ftype].longname);
141         default:
142             bPert = FALSE;
143             gmx_fatal(FARGS, "Function type %s not implemented in ip_pert",
144                       interaction_function[ftype].longname);
145     }
146
147     return bPert;
148 }
149
150 static gmx_bool ip_q_pert(int ftype, const t_iatom *ia,
151                           const t_iparams *ip, const real *qA, const real *qB)
152 {
153     /* 1-4 interactions do not have the charges stored in the iparams list,
154      * so we need a separate check for those.
155      */
156     return (ip_pert(ftype, ip+ia[0]) ||
157             (ftype == F_LJ14 && (qA[ia[1]] != qB[ia[1]] ||
158                                  qA[ia[2]] != qB[ia[2]])));
159 }
160
161 gmx_bool gmx_mtop_bondeds_free_energy(const gmx_mtop_t *mtop)
162 {
163     const gmx_ffparams_t *ffparams;
164     int                   i, ftype;
165     int                   mb;
166     t_atom               *atom;
167     t_ilist              *il;
168     t_iatom              *ia;
169     gmx_bool              bPert;
170
171     ffparams = &mtop->ffparams;
172
173     /* Loop over all the function types and compare the A/B parameters */
174     bPert = FALSE;
175     for (i = 0; i < ffparams->ntypes; i++)
176     {
177         ftype = ffparams->functype[i];
178         if (interaction_function[ftype].flags & IF_BOND)
179         {
180             if (ip_pert(ftype, &ffparams->iparams[i]))
181             {
182                 bPert = TRUE;
183             }
184         }
185     }
186
187     /* Check perturbed charges for 1-4 interactions */
188     for (mb = 0; mb < mtop->nmolblock; mb++)
189     {
190         atom = mtop->moltype[mtop->molblock[mb].type].atoms.atom;
191         il   = &mtop->moltype[mtop->molblock[mb].type].ilist[F_LJ14];
192         ia   = il->iatoms;
193         for (i = 0; i < il->nr; i += 3)
194         {
195             if (atom[ia[i+1]].q != atom[ia[i+1]].qB ||
196                 atom[ia[i+2]].q != atom[ia[i+2]].qB)
197             {
198                 bPert = TRUE;
199             }
200         }
201     }
202
203     return bPert;
204 }
205
206 void gmx_sort_ilist_fe(t_idef *idef, const real *qA, const real *qB)
207 {
208     int        ftype, nral, i, ic, ib, a;
209     t_iparams *iparams;
210     t_ilist   *ilist;
211     t_iatom   *iatoms;
212     gmx_bool   bPert;
213     t_iatom   *iabuf;
214     int        iabuf_nalloc;
215
216     if (qB == NULL)
217     {
218         qB = qA;
219     }
220
221     iabuf_nalloc = 0;
222     iabuf        = NULL;
223
224     iparams = idef->iparams;
225
226     for (ftype = 0; ftype < F_NRE; ftype++)
227     {
228         if (interaction_function[ftype].flags & IF_BOND)
229         {
230             ilist  = &idef->il[ftype];
231             iatoms = ilist->iatoms;
232             nral   = NRAL(ftype);
233             ic     = 0;
234             ib     = 0;
235             i      = 0;
236             while (i < ilist->nr)
237             {
238                 /* Check if this interaction is perturbed */
239                 if (ip_q_pert(ftype, iatoms+i, iparams, qA, qB))
240                 {
241                     /* Copy to the perturbed buffer */
242                     if (ib + 1 + nral > iabuf_nalloc)
243                     {
244                         iabuf_nalloc = over_alloc_large(ib+1+nral);
245                         srenew(iabuf, iabuf_nalloc);
246                     }
247                     for (a = 0; a < 1+nral; a++)
248                     {
249                         iabuf[ib++] = iatoms[i++];
250                     }
251                 }
252                 else
253                 {
254                     /* Copy in place */
255                     for (a = 0; a < 1+nral; a++)
256                     {
257                         iatoms[ic++] = iatoms[i++];
258                     }
259                 }
260             }
261             /* Now we now the number of non-perturbed interactions */
262             ilist->nr_nonperturbed = ic;
263
264             /* Copy the buffer with perturbed interactions to the ilist */
265             for (a = 0; a < ib; a++)
266             {
267                 iatoms[ic++] = iabuf[a];
268             }
269
270             if (debug)
271             {
272                 fprintf(debug, "%s non-pert %d pert %d\n",
273                         interaction_function[ftype].longname,
274                         ilist->nr_nonperturbed,
275                         ilist->nr-ilist->nr_nonperturbed);
276             }
277         }
278     }
279
280     sfree(iabuf);
281
282     idef->ilsort = ilsortFE_SORTED;
283 }