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)
64 bPert = (ip->harmonic.rA != ip->harmonic.rB ||
65 ip->harmonic.krA != ip->harmonic.krB);
68 bPert = (ip->morse.b0A != ip->morse.b0B ||
69 ip->morse.cbA != ip->morse.cbB ||
70 ip->morse.betaA != ip->morse.betaB);
73 bPert = (ip->restraint.lowA != ip->restraint.lowB ||
74 ip->restraint.up1A != ip->restraint.up1B ||
75 ip->restraint.up2A != ip->restraint.up2B ||
76 ip->restraint.kA != ip->restraint.kB);
79 bPert = (ip->u_b.thetaA != ip->u_b.thetaB ||
80 ip->u_b.kthetaA != ip->u_b.kthetaB ||
81 ip->u_b.r13A != ip->u_b.r13B ||
82 ip->u_b.kUBA != ip->u_b.kUBB);
88 bPert = (ip->pdihs.phiA != ip->pdihs.phiB ||
89 ip->pdihs.cpA != ip->pdihs.cpB);
93 for (i = 0; i < NR_RBDIHS; i++)
95 if (ip->rbdihs.rbcA[i] != ip->rbdihs.rbcB[i])
105 bPert = (ip->tab.kA != ip->tab.kB);
109 for (i = 0; i < DIM; i++)
111 if (ip->posres.pos0A[i] != ip->posres.pos0B[i] ||
112 ip->posres.fcA[i] != ip->posres.fcB[i])
119 bPert = ((ip->dihres.phiA != ip->dihres.phiB) ||
120 (ip->dihres.dphiA != ip->dihres.dphiB) ||
121 (ip->dihres.kfacA != ip->dihres.kfacB));
124 bPert = (ip->lj14.c6A != ip->lj14.c6B ||
125 ip->lj14.c12A != ip->lj14.c12B);
132 gmx_fatal(FARGS, "Function type %s not implemented in ip_pert",
133 interaction_function[ftype].longname);
139 static gmx_bool ip_q_pert(int ftype, const t_iatom *ia,
140 const t_iparams *ip, const real *qA, const real *qB)
142 /* 1-4 interactions do not have the charges stored in the iparams list,
143 * so we need a separate check for those.
145 return (ip_pert(ftype, ip+ia[0]) ||
146 (ftype == F_LJ14 && (qA[ia[1]] != qB[ia[1]] ||
147 qA[ia[2]] != qB[ia[2]])));
150 gmx_bool gmx_mtop_bondeds_free_energy(const gmx_mtop_t *mtop)
152 const gmx_ffparams_t *ffparams;
160 ffparams = &mtop->ffparams;
162 /* Loop over all the function types and compare the A/B parameters */
164 for (i = 0; i < ffparams->ntypes; i++)
166 ftype = ffparams->functype[i];
167 if (interaction_function[ftype].flags & IF_BOND)
169 if (ip_pert(ftype, &ffparams->iparams[i]))
176 /* Check perturbed charges for 1-4 interactions */
177 for (mb = 0; mb < mtop->nmolblock; mb++)
179 atom = mtop->moltype[mtop->molblock[mb].type].atoms.atom;
180 il = &mtop->moltype[mtop->molblock[mb].type].ilist[F_LJ14];
182 for (i = 0; i < il->nr; i += 3)
184 if (atom[ia[i+1]].q != atom[ia[i+1]].qB ||
185 atom[ia[i+2]].q != atom[ia[i+2]].qB)
195 void gmx_sort_ilist_fe(t_idef *idef, const real *qA, const real *qB)
197 int ftype, nral, i, ic, ib, a;
213 iparams = idef->iparams;
215 for (ftype = 0; ftype < F_NRE; ftype++)
217 if (interaction_function[ftype].flags & IF_BOND)
219 ilist = &idef->il[ftype];
220 iatoms = ilist->iatoms;
225 while (i < ilist->nr)
227 /* Check if this interaction is perturbed */
228 if (ip_q_pert(ftype, iatoms+i, iparams, qA, qB))
230 /* Copy to the perturbed buffer */
231 if (ib + 1 + nral > iabuf_nalloc)
233 iabuf_nalloc = over_alloc_large(ib+1+nral);
234 srenew(iabuf, iabuf_nalloc);
236 for (a = 0; a < 1+nral; a++)
238 iabuf[ib++] = iatoms[i++];
244 for (a = 0; a < 1+nral; a++)
246 iatoms[ic++] = iatoms[i++];
250 /* Now we now the number of non-perturbed interactions */
251 ilist->nr_nonperturbed = ic;
253 /* Copy the buffer with perturbed interactions to the ilist */
254 for (a = 0; a < ib; a++)
256 iatoms[ic++] = iabuf[a];
261 fprintf(debug, "%s non-pert %d pert %d\n",
262 interaction_function[ftype].longname,
263 ilist->nr_nonperturbed,
264 ilist->nr-ilist->nr_nonperturbed);
271 idef->ilsort = ilsortFE_SORTED;