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,2015,2017,2018,2019, 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.
43 #include "gromacs/topology/ifunc.h"
44 #include "gromacs/topology/topology.h"
45 #include "gromacs/utility/arrayref.h"
46 #include "gromacs/utility/fatalerror.h"
47 #include "gromacs/utility/smalloc.h"
49 static gmx_bool ip_pert(int ftype, const t_iparams* ip)
54 if (NRFPB(ftype) == 0)
67 bPert = (ip->harmonic.rA != ip->harmonic.rB || ip->harmonic.krA != ip->harmonic.krB);
70 bPert = (ip->morse.b0A != ip->morse.b0B || ip->morse.cbA != ip->morse.cbB
71 || ip->morse.betaA != ip->morse.betaB);
74 bPert = (ip->restraint.lowA != ip->restraint.lowB || 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 || ip->u_b.kthetaA != ip->u_b.kthetaB
80 || ip->u_b.r13A != ip->u_b.r13B || ip->u_b.kUBA != ip->u_b.kUBB);
86 bPert = (ip->pdihs.phiA != ip->pdihs.phiB || ip->pdihs.cpA != ip->pdihs.cpB);
90 for (i = 0; i < NR_RBDIHS; i++)
92 if (ip->rbdihs.rbcA[i] != ip->rbdihs.rbcB[i])
101 case F_TABDIHS: bPert = (ip->tab.kA != ip->tab.kB); break;
104 for (i = 0; i < DIM; i++)
106 if (ip->posres.pos0A[i] != ip->posres.pos0B[i] || ip->posres.fcA[i] != ip->posres.fcB[i])
113 bPert = ((ip->dihres.phiA != ip->dihres.phiB) || (ip->dihres.dphiA != ip->dihres.dphiB)
114 || (ip->dihres.kfacA != ip->dihres.kfacB));
117 bPert = (ip->lj14.c6A != ip->lj14.c6B || ip->lj14.c12A != ip->lj14.c12B);
119 case F_CMAP: bPert = FALSE; break;
123 gmx_fatal(FARGS, "Function type %s does not support currentely free energy calculations",
124 interaction_function[ftype].longname);
126 gmx_fatal(FARGS, "Function type %s not implemented in ip_pert",
127 interaction_function[ftype].longname);
133 static gmx_bool ip_q_pert(int ftype, const t_iatom* ia, const t_iparams* ip, const real* qA, const real* qB)
135 /* 1-4 interactions do not have the charges stored in the iparams list,
136 * so we need a separate check for those.
138 return (ip_pert(ftype, ip + ia[0])
139 || (ftype == F_LJ14 && (qA[ia[1]] != qB[ia[1]] || qA[ia[2]] != qB[ia[2]])));
142 gmx_bool gmx_mtop_bondeds_free_energy(const gmx_mtop_t* mtop)
144 const gmx_ffparams_t* ffparams = &mtop->ffparams;
146 /* Loop over all the function types and compare the A/B parameters */
147 gmx_bool bPert = FALSE;
148 for (int i = 0; i < ffparams->numTypes(); i++)
150 int ftype = ffparams->functype[i];
151 if (interaction_function[ftype].flags & IF_BOND)
153 if (ip_pert(ftype, &ffparams->iparams[i]))
160 /* Check perturbed charges for 1-4 interactions */
161 for (const gmx_molblock_t& molb : mtop->molblock)
163 const t_atom* atom = mtop->moltype[molb.type].atoms.atom;
164 const InteractionList& il = mtop->moltype[molb.type].ilist[F_LJ14];
165 gmx::ArrayRef<const int> ia = il.iatoms;
166 for (int i = 0; i < il.size(); i += 3)
168 if (atom[ia[i + 1]].q != atom[ia[i + 1]].qB || atom[ia[i + 2]].q != atom[ia[i + 2]].qB)
178 void gmx_sort_ilist_fe(t_idef* idef, const real* qA, const real* qB)
180 int ftype, nral, i, ic, ib, a;
194 const t_iparams* 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 */
235 for (a = 0; a < ib; a++)
237 iatoms[ic++] = iabuf[a];
242 fprintf(debug, "%s non-pert %d pert %d\n", interaction_function[ftype].longname,
243 ilist->nr_nonperturbed, ilist->nr - ilist->nr_nonperturbed);
250 idef->ilsort = ilsortFE_SORTED;