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