Revised wording in pdb2gmx.c, hopefully clearer now.
[gromacs/rigid-bodies.git] / src / gmxlib / topsort.c
blob128e57fb2d5b256364cf5a165c5e04c5e3123819
1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*- */
3 /*
4 *
5 * This source code is part of
6 *
7 * G R O M A C S
8 *
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
33 * And Hey:
34 * Gallium Rubidium Oxygen Manganese Argon Carbon Silicon
37 #ifdef HAVE_CONFIG_H
38 #include <config.h>
39 #endif
41 #include "typedefs.h"
42 #include "topsort.h"
43 #include "smalloc.h"
44 #include "gmx_fatal.h"
46 static gmx_bool ip_pert(int ftype,const t_iparams *ip)
48 gmx_bool bPert;
49 int i;
51 if (NRFPB(ftype) == 0)
53 return FALSE;
56 switch (ftype)
58 case F_BONDS:
59 case F_G96BONDS:
60 case F_HARMONIC:
61 case F_ANGLES:
62 case F_G96ANGLES:
63 case F_IDIHS:
64 case F_PIDIHS:
65 bPert = (ip->harmonic.rA != ip->harmonic.rB ||
66 ip->harmonic.krA != ip->harmonic.krB);
67 break;
68 case F_RESTRBONDS:
69 bPert = (ip->restraint.lowA != ip->restraint.lowB ||
70 ip->restraint.up1A != ip->restraint.up1B ||
71 ip->restraint.up2A != ip->restraint.up2B ||
72 ip->restraint.kA != ip->restraint.kB);
73 break;
74 case F_PDIHS:
75 case F_ANGRES:
76 case F_ANGRESZ:
77 bPert = (ip->pdihs.phiA != ip->pdihs.phiB ||
78 ip->pdihs.cpA != ip->pdihs.cpB);
79 break;
80 case F_RBDIHS:
81 bPert = FALSE;
82 for(i=0; i<NR_RBDIHS; i++)
84 if (ip->rbdihs.rbcA[i] != ip->rbdihs.rbcB[i])
86 bPert = TRUE;
89 break;
90 case F_TABBONDS:
91 case F_TABBONDSNC:
92 case F_TABANGLES:
93 case F_TABDIHS:
94 bPert = (ip->tab.kA != ip->tab.kB);
95 break;
96 case F_POSRES:
97 bPert = FALSE;
98 for(i=0; i<DIM; i++)
100 if (ip->posres.pos0A[i] != ip->posres.pos0B[i] ||
101 ip->posres.fcA[i] != ip->posres.fcB[i])
103 bPert = TRUE;
106 break;
107 case F_LJ14:
108 bPert = (ip->lj14.c6A != ip->lj14.c6B ||
109 ip->lj14.c12A != ip->lj14.c12B);
110 break;
111 default:
112 bPert = FALSE;
113 gmx_fatal(FARGS,"Function type %s not implemented in ip_pert",
114 interaction_function[ftype].longname);
117 return bPert;
120 static gmx_bool ip_q_pert(int ftype,const t_iatom *ia,
121 const t_iparams *ip,const real *qA,const real *qB)
123 /* 1-4 interactions do not have the charges stored in the iparams list,
124 * so we need a separate check for those.
126 return (ip_pert(ftype,ip+ia[0]) ||
127 (ftype == F_LJ14 && (qA[ia[1]] != qB[ia[1]] ||
128 qA[ia[2]] != qB[ia[2]])));
131 gmx_bool gmx_mtop_bondeds_free_energy(const gmx_mtop_t *mtop)
133 const gmx_ffparams_t *ffparams;
134 int i,ftype;
135 int mb;
136 t_atom *atom;
137 t_ilist *il;
138 t_iatom *ia;
139 gmx_bool bPert;
141 ffparams = &mtop->ffparams;
143 /* Loop over all the function types and compare the A/B parameters */
144 bPert = FALSE;
145 for(i=0; i<ffparams->ntypes; i++)
147 ftype = ffparams->functype[i];
148 if (interaction_function[ftype].flags & IF_BOND)
150 if (ip_pert(ftype,&ffparams->iparams[i]))
152 bPert = TRUE;
157 /* Check perturbed charges for 1-4 interactions */
158 for(mb=0; mb<mtop->nmolblock; mb++)
160 atom = mtop->moltype[mtop->molblock[mb].type].atoms.atom;
161 il = &mtop->moltype[mtop->molblock[mb].type].ilist[F_LJ14];
162 ia = il->iatoms;
163 for(i=0; i<il->nr; i+=3)
165 if (atom[ia[i+1]].q != atom[ia[i+1]].qB ||
166 atom[ia[i+2]].q != atom[ia[i+2]].qB)
168 bPert = TRUE;
173 return (bPert ? ilsortFE_UNSORTED : ilsortNO_FE);
176 void gmx_sort_ilist_fe(t_idef *idef,const real *qA,const real *qB)
178 int ftype,nral,i,ic,ib,a;
179 t_iparams *iparams;
180 t_ilist *ilist;
181 t_iatom *iatoms;
182 gmx_bool bPert;
183 t_iatom *iabuf;
184 int iabuf_nalloc;
186 if (qB == NULL)
188 qB = qA;
191 iabuf_nalloc = 0;
192 iabuf = NULL;
194 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;
202 nral = NRAL(ftype);
203 ic = 0;
204 ib = 0;
205 i = 0;
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++];
222 else
224 /* Copy in place */
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];
240 if (debug)
242 fprintf(debug,"%s non-pert %d pert %d\n",
243 interaction_function[ftype].longname,
244 ilist->nr_nonperturbed,
245 ilist->nr-ilist->nr_nonperturbed);
250 sfree(iabuf);
252 idef->ilsort = ilsortFE_SORTED;