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-2004, 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.
40 #include "gromacs/legacyheaders/typedefs.h"
41 #include "gromacs/math/vec.h"
42 #include "gromacs/pbcutil/pbc.h"
43 #include "gromacs/utility/smalloc.h"
44 #include "gromacs/utility/fatalerror.h"
45 #include "gromacs/legacyheaders/chargegroup.h"
48 void calc_chargegroup_radii(const gmx_mtop_t *mtop, rvec *x,
49 real *rvdw1, real *rvdw2,
50 real *rcoul1, real *rcoul2)
52 real r2v1, r2v2, r2c1, r2c2, r2;
53 int ntype, i, j, mb, m, cg, a_mol, a0, a1, a;
66 ntype = mtop->ffparams.atnr;
68 for (i = 0; i < ntype; i++)
71 for (j = 0; j < ntype; j++)
73 if (mtop->ffparams.iparams[i*ntype+j].lj.c6 != 0 ||
74 mtop->ffparams.iparams[i*ntype+j].lj.c12 != 0)
82 for (mb = 0; mb < mtop->nmolblock; mb++)
84 molb = &mtop->molblock[mb];
85 molt = &mtop->moltype[molb->type];
87 atom = molt->atoms.atom;
88 for (m = 0; m < molb->nmol; m++)
90 for (cg = 0; cg < cgs->nr; cg++)
93 a1 = cgs->index[cg+1];
97 for (a = a0; a < a1; a++)
99 rvec_inc(cen, x[a_mol+a]);
101 svmul(1.0/(a1-a0), cen, cen);
102 for (a = a0; a < a1; a++)
104 r2 = distance2(cen, x[a_mol+a]);
105 if (r2 > r2v2 && (bLJ[atom[a].type ] ||
119 (atom[a].q != 0 || atom[a].qB != 0))
134 a_mol += molb->natoms_mol;
142 *rcoul1 = sqrt(r2c1);
143 *rcoul2 = sqrt(r2c2);
146 void calc_cgcm(FILE gmx_unused *fplog, int cg0, int cg1, t_block *cgs,
147 rvec pos[], rvec cg_cm[])
149 int icg, k, k0, k1, d;
155 fprintf(fplog, "Calculating centre of geometry for charge groups %d to %d\n",
158 cgindex = cgs->index;
160 /* Compute the center of geometry for all charge groups */
161 for (icg = cg0; (icg < cg1); icg++)
168 copy_rvec(pos[k0], cg_cm[icg]);
175 for (k = k0; (k < k1); k++)
177 for (d = 0; (d < DIM); d++)
182 for (d = 0; (d < DIM); d++)
184 cg_cm[icg][d] = inv_ncg*cg[d];
190 void put_charge_groups_in_box(FILE gmx_unused *fplog, int cg0, int cg1,
191 int ePBC, matrix box, t_block *cgs,
192 rvec pos[], rvec cg_cm[])
195 int npbcdim, icg, k, k0, k1, d, e;
201 if (ePBC == epbcNONE)
203 gmx_incons("Calling put_charge_groups_in_box for a system without PBC");
207 fprintf(fplog, "Putting cgs %d to %d in box\n", cg0, cg1);
209 cgindex = cgs->index;
220 bTric = TRICLINIC(box);
222 for (icg = cg0; (icg < cg1); icg++)
224 /* First compute the center of geometry for this charge group */
231 copy_rvec(pos[k0], cg_cm[icg]);
238 for (k = k0; (k < k1); k++)
240 for (d = 0; d < DIM; d++)
245 for (d = 0; d < DIM; d++)
247 cg_cm[icg][d] = inv_ncg*cg[d];
250 /* Now check pbc for this cg */
253 for (d = npbcdim-1; d >= 0; d--)
255 while (cg_cm[icg][d] < 0)
257 for (e = d; e >= 0; e--)
259 cg_cm[icg][e] += box[d][e];
260 for (k = k0; (k < k1); k++)
262 pos[k][e] += box[d][e];
266 while (cg_cm[icg][d] >= box[d][d])
268 for (e = d; e >= 0; e--)
270 cg_cm[icg][e] -= box[d][e];
271 for (k = k0; (k < k1); k++)
273 pos[k][e] -= box[d][e];
281 for (d = 0; d < npbcdim; d++)
283 while (cg_cm[icg][d] < 0)
285 cg_cm[icg][d] += box[d][d];
286 for (k = k0; (k < k1); k++)
288 pos[k][d] += box[d][d];
291 while (cg_cm[icg][d] >= box[d][d])
293 cg_cm[icg][d] -= box[d][d];
294 for (k = k0; (k < k1); k++)
296 pos[k][d] -= box[d][d];
302 for (d = 0; (d < npbcdim); d++)
304 if ((cg_cm[icg][d] < 0) || (cg_cm[icg][d] >= box[d][d]))
306 gmx_fatal(FARGS, "cg_cm[%d] = %15f %15f %15f\n"
307 "box = %15f %15f %15f\n",
308 icg, cg_cm[icg][XX], cg_cm[icg][YY], cg_cm[icg][ZZ],
309 box[XX][XX], box[YY][YY], box[ZZ][ZZ]);