2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2008, The GROMACS development team.
6 * Copyright (c) 2013,2014, by the GROMACS development team, led by
7 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8 * and including many others, as listed in the AUTHORS file in the
9 * top-level source directory and at http://www.gromacs.org.
11 * GROMACS is free software; you can redistribute it and/or
12 * modify it under the terms of the GNU Lesser General Public License
13 * as published by the Free Software Foundation; either version 2.1
14 * of the License, or (at your option) any later version.
16 * GROMACS is distributed in the hope that it will be useful,
17 * but WITHOUT ANY WARRANTY; without even the implied warranty of
18 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 * Lesser General Public License for more details.
21 * You should have received a copy of the GNU Lesser General Public
22 * License along with GROMACS; if not, see
23 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
26 * If you want to redistribute modifications to GROMACS, please
27 * consider that scientific software is very special. Version
28 * control is crucial - bugs must be traceable. We will be happy to
29 * consider code for inclusion in the official distribution, but
30 * derived work must not be called official GROMACS. Details are found
31 * in the README & COPYING files - if they are missing, get the
32 * official version at http://www.gromacs.org.
34 * To help us fund GROMACS development, we humbly ask that you cite
35 * the research papers on the package. Check out http://www.gromacs.org.
39 #include "gromacs/topology/topsort.h"
45 #include "gromacs/legacyheaders/types/ifunc.h"
47 #include "gromacs/topology/topology.h"
48 #include "gromacs/utility/fatalerror.h"
49 #include "gromacs/utility/smalloc.h"
51 static gmx_bool ip_pert(int ftype, const t_iparams *ip)
56 if (NRFPB(ftype) == 0)
69 bPert = (ip->harmonic.rA != ip->harmonic.rB ||
70 ip->harmonic.krA != ip->harmonic.krB);
73 bPert = (ip->morse.b0A != ip->morse.b0B ||
74 ip->morse.cbA != ip->morse.cbB ||
75 ip->morse.betaA != ip->morse.betaB);
78 bPert = (ip->restraint.lowA != ip->restraint.lowB ||
79 ip->restraint.up1A != ip->restraint.up1B ||
80 ip->restraint.up2A != ip->restraint.up2B ||
81 ip->restraint.kA != ip->restraint.kB);
84 bPert = (ip->u_b.thetaA != ip->u_b.thetaB ||
85 ip->u_b.kthetaA != ip->u_b.kthetaB ||
86 ip->u_b.r13A != ip->u_b.r13B ||
87 ip->u_b.kUBA != ip->u_b.kUBB);
93 bPert = (ip->pdihs.phiA != ip->pdihs.phiB ||
94 ip->pdihs.cpA != ip->pdihs.cpB);
98 for (i = 0; i < NR_RBDIHS; i++)
100 if (ip->rbdihs.rbcA[i] != ip->rbdihs.rbcB[i])
110 bPert = (ip->tab.kA != ip->tab.kB);
114 for (i = 0; i < DIM; i++)
116 if (ip->posres.pos0A[i] != ip->posres.pos0B[i] ||
117 ip->posres.fcA[i] != ip->posres.fcB[i])
124 bPert = ((ip->dihres.phiA != ip->dihres.phiB) ||
125 (ip->dihres.dphiA != ip->dihres.dphiB) ||
126 (ip->dihres.kfacA != ip->dihres.kfacB));
129 bPert = (ip->lj14.c6A != ip->lj14.c6B ||
130 ip->lj14.c12A != ip->lj14.c12B);
139 gmx_fatal(FARGS, "Function type %s does not support currentely free energy calculations",
140 interaction_function[ftype].longname);
143 gmx_fatal(FARGS, "Function type %s not implemented in ip_pert",
144 interaction_function[ftype].longname);
150 static gmx_bool ip_q_pert(int ftype, const t_iatom *ia,
151 const t_iparams *ip, const real *qA, const real *qB)
153 /* 1-4 interactions do not have the charges stored in the iparams list,
154 * so we need a separate check for those.
156 return (ip_pert(ftype, ip+ia[0]) ||
157 (ftype == F_LJ14 && (qA[ia[1]] != qB[ia[1]] ||
158 qA[ia[2]] != qB[ia[2]])));
161 gmx_bool gmx_mtop_bondeds_free_energy(const gmx_mtop_t *mtop)
163 const gmx_ffparams_t *ffparams;
171 ffparams = &mtop->ffparams;
173 /* Loop over all the function types and compare the A/B parameters */
175 for (i = 0; i < ffparams->ntypes; i++)
177 ftype = ffparams->functype[i];
178 if (interaction_function[ftype].flags & IF_BOND)
180 if (ip_pert(ftype, &ffparams->iparams[i]))
187 /* Check perturbed charges for 1-4 interactions */
188 for (mb = 0; mb < mtop->nmolblock; mb++)
190 atom = mtop->moltype[mtop->molblock[mb].type].atoms.atom;
191 il = &mtop->moltype[mtop->molblock[mb].type].ilist[F_LJ14];
193 for (i = 0; i < il->nr; i += 3)
195 if (atom[ia[i+1]].q != atom[ia[i+1]].qB ||
196 atom[ia[i+2]].q != atom[ia[i+2]].qB)
206 void gmx_sort_ilist_fe(t_idef *idef, const real *qA, const real *qB)
208 int ftype, nral, i, ic, ib, a;
224 iparams = idef->iparams;
226 for (ftype = 0; ftype < F_NRE; ftype++)
228 if (interaction_function[ftype].flags & IF_BOND)
230 ilist = &idef->il[ftype];
231 iatoms = ilist->iatoms;
236 while (i < ilist->nr)
238 /* Check if this interaction is perturbed */
239 if (ip_q_pert(ftype, iatoms+i, iparams, qA, qB))
241 /* Copy to the perturbed buffer */
242 if (ib + 1 + nral > iabuf_nalloc)
244 iabuf_nalloc = over_alloc_large(ib+1+nral);
245 srenew(iabuf, iabuf_nalloc);
247 for (a = 0; a < 1+nral; a++)
249 iabuf[ib++] = iatoms[i++];
255 for (a = 0; a < 1+nral; a++)
257 iatoms[ic++] = iatoms[i++];
261 /* Now we now the number of non-perturbed interactions */
262 ilist->nr_nonperturbed = ic;
264 /* Copy the buffer with perturbed interactions to the ilist */
265 for (a = 0; a < ib; a++)
267 iatoms[ic++] = iabuf[a];
272 fprintf(debug, "%s non-pert %d pert %d\n",
273 interaction_function[ftype].longname,
274 ilist->nr_nonperturbed,
275 ilist->nr-ilist->nr_nonperturbed);
282 idef->ilsort = ilsortFE_SORTED;