c6c976ed5d74242bbeddb172f57beb63cd0f8d1f
[alexxy/gromacs.git] / src / gromacs / gmxlib / topsort.c
1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-  */
2
3 /*
4  *
5  *                This source code is part of
6  *
7  *                 G   R   O   M   A   C   S
8  *
9  *          GROningen MAchine for Chemical Simulations
10  *
11  * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
12  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
13  * Copyright (c) 2001-2008, The GROMACS development team,
14  * check out http://www.gromacs.org for more information.
15
16  * This program is free software; you can redistribute it and/or
17  * modify it under the terms of the GNU General Public License
18  * as published by the Free Software Foundation; either version 2
19  * of the License, or (at your option) any later version.
20  *
21  * If you want to redistribute modifications, please consider that
22  * scientific software is very special. Version control is crucial -
23  * bugs must be traceable. We will be happy to consider code for
24  * inclusion in the official distribution, but derived work must not
25  * be called official GROMACS. Details are found in the README & COPYING
26  * files - if they are missing, get the official version at www.gromacs.org.
27  *
28  * To help us fund GROMACS development, we humbly ask that you cite
29  * the papers on the package - you can find them in the top README file.
30  *
31  * For more info, check our website at http://www.gromacs.org
32  *
33  * And Hey:
34  * Gallium Rubidium Oxygen Manganese Argon Carbon Silicon
35  */
36
37 #ifdef HAVE_CONFIG_H
38 #include <config.h>
39 #endif
40
41 #include "typedefs.h"
42 #include "topsort.h"
43 #include "smalloc.h"
44 #include "gmx_fatal.h"
45
46 static gmx_bool ip_pert(int ftype, const t_iparams *ip)
47 {
48     gmx_bool bPert;
49     int      i;
50
51     if (NRFPB(ftype) == 0)
52     {
53         return FALSE;
54     }
55
56     switch (ftype)
57     {
58         case F_BONDS:
59         case F_G96BONDS:
60         case F_HARMONIC:
61         case F_ANGLES:
62         case F_G96ANGLES:
63         case F_IDIHS:
64             bPert = (ip->harmonic.rA  != ip->harmonic.rB ||
65                      ip->harmonic.krA != ip->harmonic.krB);
66             break;
67         case F_MORSE:
68             bPert = (ip->morse.b0A  != ip->morse.b0B ||
69                      ip->morse.cbA  != ip->morse.cbB ||
70                      ip->morse.betaA  != ip->morse.betaB);
71             break;
72         case F_RESTRBONDS:
73             bPert = (ip->restraint.lowA  != ip->restraint.lowB ||
74                      ip->restraint.up1A  != ip->restraint.up1B ||
75                      ip->restraint.up2A  != ip->restraint.up2B ||
76                      ip->restraint.kA    != ip->restraint.kB);
77             break;
78         case F_UREY_BRADLEY:
79             bPert = (ip->u_b.thetaA  != ip->u_b.thetaB  ||
80                      ip->u_b.kthetaA != ip->u_b.kthetaB ||
81                      ip->u_b.r13A    != ip->u_b.r13B    ||
82                      ip->u_b.kUBA    != ip->u_b.kUBB);
83             break;
84         case F_PDIHS:
85         case F_PIDIHS:
86         case F_ANGRES:
87         case F_ANGRESZ:
88             bPert = (ip->pdihs.phiA != ip->pdihs.phiB ||
89                      ip->pdihs.cpA  != ip->pdihs.cpB);
90             break;
91         case F_RBDIHS:
92             bPert = FALSE;
93             for (i = 0; i < NR_RBDIHS; i++)
94             {
95                 if (ip->rbdihs.rbcA[i] != ip->rbdihs.rbcB[i])
96                 {
97                     bPert = TRUE;
98                 }
99             }
100             break;
101         case F_TABBONDS:
102         case F_TABBONDSNC:
103         case F_TABANGLES:
104         case F_TABDIHS:
105             bPert = (ip->tab.kA != ip->tab.kB);
106             break;
107         case F_POSRES:
108             bPert = FALSE;
109             for (i = 0; i < DIM; i++)
110             {
111                 if (ip->posres.pos0A[i] != ip->posres.pos0B[i] ||
112                     ip->posres.fcA[i]   != ip->posres.fcB[i])
113                 {
114                     bPert = TRUE;
115                 }
116             }
117             break;
118         case F_DIHRES:
119             bPert = ((ip->dihres.phiA != ip->dihres.phiB) ||
120                      (ip->dihres.dphiA != ip->dihres.dphiB) ||
121                      (ip->dihres.kfacA != ip->dihres.kfacB));
122             break;
123         case F_LJ14:
124             bPert = (ip->lj14.c6A  != ip->lj14.c6B ||
125                      ip->lj14.c12A != ip->lj14.c12B);
126             break;
127         case F_CMAP:
128             bPert = FALSE;
129             break;
130         default:
131             bPert = FALSE;
132             gmx_fatal(FARGS, "Function type %s not implemented in ip_pert",
133                       interaction_function[ftype].longname);
134     }
135
136     return bPert;
137 }
138
139 static gmx_bool ip_q_pert(int ftype, const t_iatom *ia,
140                           const t_iparams *ip, const real *qA, const real *qB)
141 {
142     /* 1-4 interactions do not have the charges stored in the iparams list,
143      * so we need a separate check for those.
144      */
145     return (ip_pert(ftype, ip+ia[0]) ||
146             (ftype == F_LJ14 && (qA[ia[1]] != qB[ia[1]] ||
147                                  qA[ia[2]] != qB[ia[2]])));
148 }
149
150 gmx_bool gmx_mtop_bondeds_free_energy(const gmx_mtop_t *mtop)
151 {
152     const gmx_ffparams_t *ffparams;
153     int                   i, ftype;
154     int                   mb;
155     t_atom               *atom;
156     t_ilist              *il;
157     t_iatom              *ia;
158     gmx_bool              bPert;
159
160     ffparams = &mtop->ffparams;
161
162     /* Loop over all the function types and compare the A/B parameters */
163     bPert = FALSE;
164     for (i = 0; i < ffparams->ntypes; i++)
165     {
166         ftype = ffparams->functype[i];
167         if (interaction_function[ftype].flags & IF_BOND)
168         {
169             if (ip_pert(ftype, &ffparams->iparams[i]))
170             {
171                 bPert = TRUE;
172             }
173         }
174     }
175
176     /* Check perturbed charges for 1-4 interactions */
177     for (mb = 0; mb < mtop->nmolblock; mb++)
178     {
179         atom = mtop->moltype[mtop->molblock[mb].type].atoms.atom;
180         il   = &mtop->moltype[mtop->molblock[mb].type].ilist[F_LJ14];
181         ia   = il->iatoms;
182         for (i = 0; i < il->nr; i += 3)
183         {
184             if (atom[ia[i+1]].q != atom[ia[i+1]].qB ||
185                 atom[ia[i+2]].q != atom[ia[i+2]].qB)
186             {
187                 bPert = TRUE;
188             }
189         }
190     }
191
192     return bPert;
193 }
194
195 void gmx_sort_ilist_fe(t_idef *idef, const real *qA, const real *qB)
196 {
197     int        ftype, nral, i, ic, ib, a;
198     t_iparams *iparams;
199     t_ilist   *ilist;
200     t_iatom   *iatoms;
201     gmx_bool   bPert;
202     t_iatom   *iabuf;
203     int        iabuf_nalloc;
204
205     if (qB == NULL)
206     {
207         qB = qA;
208     }
209
210     iabuf_nalloc = 0;
211     iabuf        = NULL;
212
213     iparams = idef->iparams;
214
215     for (ftype = 0; ftype < F_NRE; ftype++)
216     {
217         if (interaction_function[ftype].flags & IF_BOND)
218         {
219             ilist  = &idef->il[ftype];
220             iatoms = ilist->iatoms;
221             nral   = NRAL(ftype);
222             ic     = 0;
223             ib     = 0;
224             i      = 0;
225             while (i < ilist->nr)
226             {
227                 /* Check if this interaction is perturbed */
228                 if (ip_q_pert(ftype, iatoms+i, iparams, qA, qB))
229                 {
230                     /* Copy to the perturbed buffer */
231                     if (ib + 1 + nral > iabuf_nalloc)
232                     {
233                         iabuf_nalloc = over_alloc_large(ib+1+nral);
234                         srenew(iabuf, iabuf_nalloc);
235                     }
236                     for (a = 0; a < 1+nral; a++)
237                     {
238                         iabuf[ib++] = iatoms[i++];
239                     }
240                 }
241                 else
242                 {
243                     /* Copy in place */
244                     for (a = 0; a < 1+nral; a++)
245                     {
246                         iatoms[ic++] = iatoms[i++];
247                     }
248                 }
249             }
250             /* Now we now the number of non-perturbed interactions */
251             ilist->nr_nonperturbed = ic;
252
253             /* Copy the buffer with perturbed interactions to the ilist */
254             for (a = 0; a < ib; a++)
255             {
256                 iatoms[ic++] = iabuf[a];
257             }
258
259             if (debug)
260             {
261                 fprintf(debug, "%s non-pert %d pert %d\n",
262                         interaction_function[ftype].longname,
263                         ilist->nr_nonperturbed,
264                         ilist->nr-ilist->nr_nonperturbed);
265             }
266         }
267     }
268
269     sfree(iabuf);
270
271     idef->ilsort = ilsortFE_SORTED;
272 }