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