Merge branch 'release-4-6'
[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_PDIHS:
79     case F_PIDIHS:
80     case F_ANGRES:
81     case F_ANGRESZ:
82         bPert = (ip->pdihs.phiA != ip->pdihs.phiB ||
83                  ip->pdihs.cpA  != ip->pdihs.cpB);
84         break;
85     case F_RBDIHS:
86         bPert = FALSE;
87         for(i=0; i<NR_RBDIHS; i++)
88         {
89             if (ip->rbdihs.rbcA[i] != ip->rbdihs.rbcB[i])
90             {
91                 bPert = TRUE;
92             }
93         }
94         break;
95     case F_TABBONDS:
96     case F_TABBONDSNC:
97     case F_TABANGLES:
98     case F_TABDIHS:
99         bPert = (ip->tab.kA != ip->tab.kB);
100         break;
101     case F_POSRES:
102         bPert = FALSE;
103         for(i=0; i<DIM; i++)
104         {
105             if (ip->posres.pos0A[i] != ip->posres.pos0B[i] ||
106                 ip->posres.fcA[i]   != ip->posres.fcB[i])
107             {
108                 bPert = TRUE;
109             }
110         }
111         break;
112     case F_DIHRES:
113         bPert = ((ip->dihres.phiA != ip->dihres.phiB) ||
114                  (ip->dihres.dphiA != ip->dihres.dphiB) ||
115                  (ip->dihres.kfacA != ip->dihres.kfacB));
116         break;
117    case F_LJ14:
118         bPert = (ip->lj14.c6A  != ip->lj14.c6B ||
119                  ip->lj14.c12A != ip->lj14.c12B);
120         break;
121     case F_CMAP:
122         bPert = FALSE;
123         break;
124     default:
125         bPert = FALSE;
126         gmx_fatal(FARGS,"Function type %s not implemented in ip_pert",
127                   interaction_function[ftype].longname);
128     }
129
130     return bPert;
131 }
132
133 static gmx_bool ip_q_pert(int ftype,const t_iatom *ia,
134                       const t_iparams *ip,const real *qA,const real *qB)
135 {
136     /* 1-4 interactions do not have the charges stored in the iparams list,
137      * so we need a separate check for those.
138      */
139     return (ip_pert(ftype,ip+ia[0]) || 
140             (ftype == F_LJ14 && (qA[ia[1]] != qB[ia[1]] ||
141                                  qA[ia[2]] != qB[ia[2]])));
142 }
143
144 gmx_bool gmx_mtop_bondeds_free_energy(const gmx_mtop_t *mtop)
145 {
146     const gmx_ffparams_t *ffparams;
147     int  i,ftype;
148     int  mb;
149     t_atom  *atom;
150     t_ilist *il;
151     t_iatom *ia;
152     gmx_bool bPert;
153
154     ffparams = &mtop->ffparams;
155     
156     /* Loop over all the function types and compare the A/B parameters */
157     bPert = FALSE;
158     for(i=0; i<ffparams->ntypes; i++)
159     {
160         ftype = ffparams->functype[i];
161         if (interaction_function[ftype].flags & IF_BOND)
162         {
163             if (ip_pert(ftype,&ffparams->iparams[i]))
164             {
165                 bPert = TRUE;
166             }
167         }
168     }
169
170     /* Check perturbed charges for 1-4 interactions */
171     for(mb=0; mb<mtop->nmolblock; mb++)
172     {
173         atom = mtop->moltype[mtop->molblock[mb].type].atoms.atom;
174         il   = &mtop->moltype[mtop->molblock[mb].type].ilist[F_LJ14];
175         ia   = il->iatoms;
176         for(i=0; i<il->nr; i+=3)
177         {
178             if (atom[ia[i+1]].q != atom[ia[i+1]].qB ||
179                 atom[ia[i+2]].q != atom[ia[i+2]].qB)
180             {
181                 bPert = TRUE;
182             }
183         }
184     }
185
186     return (bPert ? ilsortFE_UNSORTED : ilsortNO_FE);
187 }
188
189 void gmx_sort_ilist_fe(t_idef *idef,const real *qA,const real *qB)
190 {
191     int  ftype,nral,i,ic,ib,a;
192     t_iparams *iparams;
193     t_ilist *ilist;
194     t_iatom *iatoms;
195     gmx_bool bPert;
196     t_iatom *iabuf;
197     int  iabuf_nalloc;
198
199     if (qB == NULL)
200     {
201         qB = qA;
202     }
203
204     iabuf_nalloc = 0;
205     iabuf        = NULL;
206     
207     iparams = idef->iparams;
208
209     for(ftype=0; ftype<F_NRE; ftype++)
210     {
211         if (interaction_function[ftype].flags & IF_BOND)
212         {
213             ilist = &idef->il[ftype];
214             iatoms = ilist->iatoms;
215             nral  = NRAL(ftype);
216             ic = 0;
217             ib = 0;
218             i  = 0;
219             while (i < ilist->nr)
220             {
221                 /* Check if this interaction is perturbed */
222                 if (ip_q_pert(ftype,iatoms+i,iparams,qA,qB))
223                 {
224                     /* Copy to the perturbed buffer */
225                     if (ib + 1 + nral > iabuf_nalloc)
226                     {
227                         iabuf_nalloc = over_alloc_large(ib+1+nral);
228                         srenew(iabuf,iabuf_nalloc);
229                     }
230                     for(a=0; a<1+nral; a++)
231                     {
232                         iabuf[ib++] = iatoms[i++];
233                     }
234                 }
235                 else
236                 {
237                     /* Copy in place */
238                     for(a=0; a<1+nral; a++)
239                     {
240                         iatoms[ic++] = iatoms[i++];
241                     }
242                 }
243             }
244             /* Now we now the number of non-perturbed interactions */
245             ilist->nr_nonperturbed = ic;
246             
247             /* Copy the buffer with perturbed interactions to the ilist */
248             for(a=0; a<ib; a++)
249             {
250                 iatoms[ic++] = iabuf[a];
251             }
252
253             if (debug)
254             {
255                 fprintf(debug,"%s non-pert %d pert %d\n",
256                         interaction_function[ftype].longname,
257                         ilist->nr_nonperturbed,
258                         ilist->nr-ilist->nr_nonperturbed);
259             }
260         }
261     }
262
263     sfree(iabuf);
264
265     idef->ilsort = ilsortFE_SORTED;
266 }