Remove all unnecessary HAVE_CONFIG_H
[gromacs.git] / src / gromacs / topology / topsort.c
blobab2669c6781da49c50bfa1c99148b610f411c613
1 /*
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.
37 #include "gromacs/topology/topsort.h"
39 #include "config.h"
41 #include <stdio.h>
43 #include "gromacs/legacyheaders/types/ifunc.h"
45 #include "gromacs/topology/topology.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)
51 gmx_bool bPert;
52 int i;
54 if (NRFPB(ftype) == 0)
56 return FALSE;
59 switch (ftype)
61 case F_BONDS:
62 case F_G96BONDS:
63 case F_HARMONIC:
64 case F_ANGLES:
65 case F_G96ANGLES:
66 case F_IDIHS:
67 bPert = (ip->harmonic.rA != ip->harmonic.rB ||
68 ip->harmonic.krA != ip->harmonic.krB);
69 break;
70 case F_MORSE:
71 bPert = (ip->morse.b0A != ip->morse.b0B ||
72 ip->morse.cbA != ip->morse.cbB ||
73 ip->morse.betaA != ip->morse.betaB);
74 break;
75 case F_RESTRBONDS:
76 bPert = (ip->restraint.lowA != ip->restraint.lowB ||
77 ip->restraint.up1A != ip->restraint.up1B ||
78 ip->restraint.up2A != ip->restraint.up2B ||
79 ip->restraint.kA != ip->restraint.kB);
80 break;
81 case F_UREY_BRADLEY:
82 bPert = (ip->u_b.thetaA != ip->u_b.thetaB ||
83 ip->u_b.kthetaA != ip->u_b.kthetaB ||
84 ip->u_b.r13A != ip->u_b.r13B ||
85 ip->u_b.kUBA != ip->u_b.kUBB);
86 break;
87 case F_PDIHS:
88 case F_PIDIHS:
89 case F_ANGRES:
90 case F_ANGRESZ:
91 bPert = (ip->pdihs.phiA != ip->pdihs.phiB ||
92 ip->pdihs.cpA != ip->pdihs.cpB);
93 break;
94 case F_RBDIHS:
95 bPert = FALSE;
96 for (i = 0; i < NR_RBDIHS; i++)
98 if (ip->rbdihs.rbcA[i] != ip->rbdihs.rbcB[i])
100 bPert = TRUE;
103 break;
104 case F_TABBONDS:
105 case F_TABBONDSNC:
106 case F_TABANGLES:
107 case F_TABDIHS:
108 bPert = (ip->tab.kA != ip->tab.kB);
109 break;
110 case F_POSRES:
111 bPert = FALSE;
112 for (i = 0; i < DIM; i++)
114 if (ip->posres.pos0A[i] != ip->posres.pos0B[i] ||
115 ip->posres.fcA[i] != ip->posres.fcB[i])
117 bPert = TRUE;
120 break;
121 case F_DIHRES:
122 bPert = ((ip->dihres.phiA != ip->dihres.phiB) ||
123 (ip->dihres.dphiA != ip->dihres.dphiB) ||
124 (ip->dihres.kfacA != ip->dihres.kfacB));
125 break;
126 case F_LJ14:
127 bPert = (ip->lj14.c6A != ip->lj14.c6B ||
128 ip->lj14.c12A != ip->lj14.c12B);
129 break;
130 case F_CMAP:
131 bPert = FALSE;
132 break;
133 case F_RESTRANGLES:
134 case F_RESTRDIHS:
135 case F_CBTDIHS:
136 bPert = FALSE;
137 gmx_fatal(FARGS, "Function type %s does not support currentely free energy calculations",
138 interaction_function[ftype].longname);
139 default:
140 bPert = FALSE;
141 gmx_fatal(FARGS, "Function type %s not implemented in ip_pert",
142 interaction_function[ftype].longname);
145 return bPert;
148 static gmx_bool ip_q_pert(int ftype, const t_iatom *ia,
149 const t_iparams *ip, const real *qA, const real *qB)
151 /* 1-4 interactions do not have the charges stored in the iparams list,
152 * so we need a separate check for those.
154 return (ip_pert(ftype, ip+ia[0]) ||
155 (ftype == F_LJ14 && (qA[ia[1]] != qB[ia[1]] ||
156 qA[ia[2]] != qB[ia[2]])));
159 gmx_bool gmx_mtop_bondeds_free_energy(const gmx_mtop_t *mtop)
161 const gmx_ffparams_t *ffparams;
162 int i, ftype;
163 int mb;
164 t_atom *atom;
165 t_ilist *il;
166 t_iatom *ia;
167 gmx_bool bPert;
169 ffparams = &mtop->ffparams;
171 /* Loop over all the function types and compare the A/B parameters */
172 bPert = FALSE;
173 for (i = 0; i < ffparams->ntypes; i++)
175 ftype = ffparams->functype[i];
176 if (interaction_function[ftype].flags & IF_BOND)
178 if (ip_pert(ftype, &ffparams->iparams[i]))
180 bPert = TRUE;
185 /* Check perturbed charges for 1-4 interactions */
186 for (mb = 0; mb < mtop->nmolblock; mb++)
188 atom = mtop->moltype[mtop->molblock[mb].type].atoms.atom;
189 il = &mtop->moltype[mtop->molblock[mb].type].ilist[F_LJ14];
190 ia = il->iatoms;
191 for (i = 0; i < il->nr; i += 3)
193 if (atom[ia[i+1]].q != atom[ia[i+1]].qB ||
194 atom[ia[i+2]].q != atom[ia[i+2]].qB)
196 bPert = TRUE;
201 return bPert;
204 void gmx_sort_ilist_fe(t_idef *idef, const real *qA, const real *qB)
206 int ftype, nral, i, ic, ib, a;
207 t_iparams *iparams;
208 t_ilist *ilist;
209 t_iatom *iatoms;
210 gmx_bool bPert;
211 t_iatom *iabuf;
212 int iabuf_nalloc;
214 if (qB == NULL)
216 qB = qA;
219 iabuf_nalloc = 0;
220 iabuf = NULL;
222 iparams = idef->iparams;
224 for (ftype = 0; ftype < F_NRE; ftype++)
226 if (interaction_function[ftype].flags & IF_BOND)
228 ilist = &idef->il[ftype];
229 iatoms = ilist->iatoms;
230 nral = NRAL(ftype);
231 ic = 0;
232 ib = 0;
233 i = 0;
234 while (i < ilist->nr)
236 /* Check if this interaction is perturbed */
237 if (ip_q_pert(ftype, iatoms+i, iparams, qA, qB))
239 /* Copy to the perturbed buffer */
240 if (ib + 1 + nral > iabuf_nalloc)
242 iabuf_nalloc = over_alloc_large(ib+1+nral);
243 srenew(iabuf, iabuf_nalloc);
245 for (a = 0; a < 1+nral; a++)
247 iabuf[ib++] = iatoms[i++];
250 else
252 /* Copy in place */
253 for (a = 0; a < 1+nral; a++)
255 iatoms[ic++] = iatoms[i++];
259 /* Now we now the number of non-perturbed interactions */
260 ilist->nr_nonperturbed = ic;
262 /* Copy the buffer with perturbed interactions to the ilist */
263 for (a = 0; a < ib; a++)
265 iatoms[ic++] = iabuf[a];
268 if (debug)
270 fprintf(debug, "%s non-pert %d pert %d\n",
271 interaction_function[ftype].longname,
272 ilist->nr_nonperturbed,
273 ilist->nr-ilist->nr_nonperturbed);
278 sfree(iabuf);
280 idef->ilsort = ilsortFE_SORTED;