1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*- */
5 * This source code is part of
9 * GROningen MAchine for Chemical Simulations
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.
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.
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.
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.
31 * For more info, check our website at http://www.gromacs.org
34 * Gallium Rubidium Oxygen Manganese Argon Carbon Silicon
44 #include "gmx_fatal.h"
46 static gmx_bool
ip_pert(int ftype
,const t_iparams
*ip
)
51 if (NRFPB(ftype
) == 0)
65 bPert
= (ip
->harmonic
.rA
!= ip
->harmonic
.rB
||
66 ip
->harmonic
.krA
!= ip
->harmonic
.krB
);
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
);
77 bPert
= (ip
->pdihs
.phiA
!= ip
->pdihs
.phiB
||
78 ip
->pdihs
.cpA
!= ip
->pdihs
.cpB
);
82 for(i
=0; i
<NR_RBDIHS
; i
++)
84 if (ip
->rbdihs
.rbcA
[i
] != ip
->rbdihs
.rbcB
[i
])
94 bPert
= (ip
->tab
.kA
!= ip
->tab
.kB
);
100 if (ip
->posres
.pos0A
[i
] != ip
->posres
.pos0B
[i
] ||
101 ip
->posres
.fcA
[i
] != ip
->posres
.fcB
[i
])
108 bPert
= (ip
->lj14
.c6A
!= ip
->lj14
.c6B
||
109 ip
->lj14
.c12A
!= ip
->lj14
.c12B
);
113 gmx_fatal(FARGS
,"Function type %s not implemented in ip_pert",
114 interaction_function
[ftype
].longname
);
120 static gmx_bool
ip_q_pert(int ftype
,const t_iatom
*ia
,
121 const t_iparams
*ip
,const real
*qA
,const real
*qB
)
123 /* 1-4 interactions do not have the charges stored in the iparams list,
124 * so we need a separate check for those.
126 return (ip_pert(ftype
,ip
+ia
[0]) ||
127 (ftype
== F_LJ14
&& (qA
[ia
[1]] != qB
[ia
[1]] ||
128 qA
[ia
[2]] != qB
[ia
[2]])));
131 gmx_bool
gmx_mtop_bondeds_free_energy(const gmx_mtop_t
*mtop
)
133 const gmx_ffparams_t
*ffparams
;
141 ffparams
= &mtop
->ffparams
;
143 /* Loop over all the function types and compare the A/B parameters */
145 for(i
=0; i
<ffparams
->ntypes
; i
++)
147 ftype
= ffparams
->functype
[i
];
148 if (interaction_function
[ftype
].flags
& IF_BOND
)
150 if (ip_pert(ftype
,&ffparams
->iparams
[i
]))
157 /* Check perturbed charges for 1-4 interactions */
158 for(mb
=0; mb
<mtop
->nmolblock
; mb
++)
160 atom
= mtop
->moltype
[mtop
->molblock
[mb
].type
].atoms
.atom
;
161 il
= &mtop
->moltype
[mtop
->molblock
[mb
].type
].ilist
[F_LJ14
];
163 for(i
=0; i
<il
->nr
; i
+=3)
165 if (atom
[ia
[i
+1]].q
!= atom
[ia
[i
+1]].qB
||
166 atom
[ia
[i
+2]].q
!= atom
[ia
[i
+2]].qB
)
173 return (bPert
? ilsortFE_UNSORTED
: ilsortNO_FE
);
176 void gmx_sort_ilist_fe(t_idef
*idef
,const real
*qA
,const real
*qB
)
178 int ftype
,nral
,i
,ic
,ib
,a
;
194 iparams
= idef
->iparams
;
196 for(ftype
=0; ftype
<F_NRE
; ftype
++)
198 if (interaction_function
[ftype
].flags
& IF_BOND
)
200 ilist
= &idef
->il
[ftype
];
201 iatoms
= ilist
->iatoms
;
206 while (i
< ilist
->nr
)
208 /* Check if this interaction is perturbed */
209 if (ip_q_pert(ftype
,iatoms
+i
,iparams
,qA
,qB
))
211 /* Copy to the perturbed buffer */
212 if (ib
+ 1 + nral
> iabuf_nalloc
)
214 iabuf_nalloc
= over_alloc_large(ib
+1+nral
);
215 srenew(iabuf
,iabuf_nalloc
);
217 for(a
=0; a
<1+nral
; a
++)
219 iabuf
[ib
++] = iatoms
[i
++];
225 for(a
=0; a
<1+nral
; a
++)
227 iatoms
[ic
++] = iatoms
[i
++];
231 /* Now we now the number of non-perturbed interactions */
232 ilist
->nr_nonperturbed
= ic
;
234 /* Copy the buffer with perturbed interactions to the ilist */
237 iatoms
[ic
++] = iabuf
[a
];
242 fprintf(debug
,"%s non-pert %d pert %d\n",
243 interaction_function
[ftype
].longname
,
244 ilist
->nr_nonperturbed
,
245 ilist
->nr
-ilist
->nr_nonperturbed
);
252 idef
->ilsort
= ilsortFE_SORTED
;