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