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.
46 #include "gromacs/utility/smalloc.h"
47 #include "gmx_fatal.h"
48 #include "chargegroup.h"
51 void calc_chargegroup_radii(const gmx_mtop_t *mtop, rvec *x,
52 real *rvdw1, real *rvdw2,
53 real *rcoul1, real *rcoul2)
55 real r2v1, r2v2, r2c1, r2c2, r2;
56 int ntype, i, j, mb, m, cg, a_mol, a0, a1, a;
69 ntype = mtop->ffparams.atnr;
71 for (i = 0; i < ntype; i++)
74 for (j = 0; j < ntype; j++)
76 if (mtop->ffparams.iparams[i*ntype+j].lj.c6 != 0 ||
77 mtop->ffparams.iparams[i*ntype+j].lj.c12 != 0)
85 for (mb = 0; mb < mtop->nmolblock; mb++)
87 molb = &mtop->molblock[mb];
88 molt = &mtop->moltype[molb->type];
90 atom = molt->atoms.atom;
91 for (m = 0; m < molb->nmol; m++)
93 for (cg = 0; cg < cgs->nr; cg++)
96 a1 = cgs->index[cg+1];
100 for (a = a0; a < a1; a++)
102 rvec_inc(cen, x[a_mol+a]);
104 svmul(1.0/(a1-a0), cen, cen);
105 for (a = a0; a < a1; a++)
107 r2 = distance2(cen, x[a_mol+a]);
108 if (r2 > r2v2 && (bLJ[atom[a].type ] ||
122 (atom[a].q != 0 || atom[a].qB != 0))
137 a_mol += molb->natoms_mol;
145 *rcoul1 = sqrt(r2c1);
146 *rcoul2 = sqrt(r2c2);
149 void calc_cgcm(FILE gmx_unused *fplog, int cg0, int cg1, t_block *cgs,
150 rvec pos[], rvec cg_cm[])
152 int icg, k, k0, k1, d;
158 fprintf(fplog, "Calculating centre of geometry for charge groups %d to %d\n",
161 cgindex = cgs->index;
163 /* Compute the center of geometry for all charge groups */
164 for (icg = cg0; (icg < cg1); icg++)
171 copy_rvec(pos[k0], cg_cm[icg]);
178 for (k = k0; (k < k1); k++)
180 for (d = 0; (d < DIM); d++)
185 for (d = 0; (d < DIM); d++)
187 cg_cm[icg][d] = inv_ncg*cg[d];
193 void put_charge_groups_in_box(FILE gmx_unused *fplog, int cg0, int cg1,
194 int ePBC, matrix box, t_block *cgs,
195 rvec pos[], rvec cg_cm[])
198 int npbcdim, icg, k, k0, k1, d, e;
204 if (ePBC == epbcNONE)
206 gmx_incons("Calling put_charge_groups_in_box for a system without PBC");
210 fprintf(fplog, "Putting cgs %d to %d in box\n", cg0, cg1);
212 cgindex = cgs->index;
223 bTric = TRICLINIC(box);
225 for (icg = cg0; (icg < cg1); icg++)
227 /* First compute the center of geometry for this charge group */
234 copy_rvec(pos[k0], cg_cm[icg]);
241 for (k = k0; (k < k1); k++)
243 for (d = 0; d < DIM; d++)
248 for (d = 0; d < DIM; d++)
250 cg_cm[icg][d] = inv_ncg*cg[d];
253 /* Now check pbc for this cg */
256 for (d = npbcdim-1; d >= 0; d--)
258 while (cg_cm[icg][d] < 0)
260 for (e = d; e >= 0; e--)
262 cg_cm[icg][e] += box[d][e];
263 for (k = k0; (k < k1); k++)
265 pos[k][e] += box[d][e];
269 while (cg_cm[icg][d] >= box[d][d])
271 for (e = d; e >= 0; e--)
273 cg_cm[icg][e] -= box[d][e];
274 for (k = k0; (k < k1); k++)
276 pos[k][e] -= box[d][e];
284 for (d = 0; d < npbcdim; d++)
286 while (cg_cm[icg][d] < 0)
288 cg_cm[icg][d] += box[d][d];
289 for (k = k0; (k < k1); k++)
291 pos[k][d] += box[d][d];
294 while (cg_cm[icg][d] >= box[d][d])
296 cg_cm[icg][d] -= box[d][d];
297 for (k = k0; (k < k1); k++)
299 pos[k][d] -= box[d][d];
305 for (d = 0; (d < npbcdim); d++)
307 if ((cg_cm[icg][d] < 0) || (cg_cm[icg][d] >= box[d][d]))
309 gmx_fatal(FARGS, "cg_cm[%d] = %15f %15f %15f\n"
310 "box = %15f %15f %15f\n",
311 icg, cg_cm[icg][XX], cg_cm[icg][YY], cg_cm[icg][ZZ],
312 box[XX][XX], box[YY][YY], box[ZZ][ZZ]);