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->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);
76 bPert = (ip->pdihs.phiA != ip->pdihs.phiB ||
77 ip->pdihs.cpA != ip->pdihs.cpB);
81 for(i=0; i<NR_RBDIHS; i++)
83 if (ip->rbdihs.rbcA[i] != ip->rbdihs.rbcB[i])
93 bPert = (ip->tab.kA != ip->tab.kB);
99 if (ip->posres.pos0A[i] != ip->posres.pos0B[i] ||
100 ip->posres.fcA[i] != ip->posres.fcB[i])
107 bPert = (ip->lj14.c6A != ip->lj14.c6B ||
108 ip->lj14.c12A != ip->lj14.c12B);
112 gmx_fatal(FARGS,"Function type %s not implemented in ip_pert",
113 interaction_function[ftype].longname);
119 static gmx_bool ip_q_pert(int ftype,const t_iatom *ia,
120 const t_iparams *ip,const real *qA,const real *qB)
122 /* 1-4 interactions do not have the charges stored in the iparams list,
123 * so we need a separate check for those.
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]])));
130 gmx_bool gmx_mtop_bondeds_free_energy(const gmx_mtop_t *mtop)
132 const gmx_ffparams_t *ffparams;
140 ffparams = &mtop->ffparams;
142 /* Loop over all the function types and compare the A/B parameters */
144 for(i=0; i<ffparams->ntypes; i++)
146 ftype = ffparams->functype[i];
147 if (interaction_function[ftype].flags & IF_BOND)
149 if (ip_pert(ftype,&ffparams->iparams[i]))
156 /* Check perturbed charges for 1-4 interactions */
157 for(mb=0; mb<mtop->nmolblock; mb++)
159 atom = mtop->moltype[mtop->molblock[mb].type].atoms.atom;
160 il = &mtop->moltype[mtop->molblock[mb].type].ilist[F_LJ14];
162 for(i=0; i<il->nr; i+=3)
164 if (atom[ia[i+1]].q != atom[ia[i+1]].qB ||
165 atom[ia[i+2]].q != atom[ia[i+2]].qB)
172 return (bPert ? ilsortFE_UNSORTED : ilsortNO_FE);
175 void gmx_sort_ilist_fe(t_idef *idef,const real *qA,const real *qB)
177 int ftype,nral,i,ic,ib,a;
193 iparams = idef->iparams;
195 for(ftype=0; ftype<F_NRE; ftype++)
197 if (interaction_function[ftype].flags & IF_BOND)
199 ilist = &idef->il[ftype];
200 iatoms = ilist->iatoms;
205 while (i < ilist->nr)
207 /* Check if this interaction is perturbed */
208 if (ip_q_pert(ftype,iatoms+i,iparams,qA,qB))
210 /* Copy to the perturbed buffer */
211 if (ib + 1 + nral > iabuf_nalloc)
213 iabuf_nalloc = over_alloc_large(ib+1+nral);
214 srenew(iabuf,iabuf_nalloc);
216 for(a=0; a<1+nral; a++)
218 iabuf[ib++] = iatoms[i++];
224 for(a=0; a<1+nral; a++)
226 iatoms[ic++] = iatoms[i++];
230 /* Now we now the number of non-perturbed interactions */
231 ilist->nr_nonperturbed = ic;
233 /* Copy the buffer with perturbed interactions to the ilist */
236 iatoms[ic++] = iabuf[a];
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);
251 idef->ilsort = ilsortFE_SORTED;